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 if (file_version >= 75)
277 gmx_fio_do_int(fio,rotg->PotAngle_nstep);
278 gmx_fio_do_real(fio,rotg->PotAngle_step);
282 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
286 gmx_fio_do_int(fio,rot->ngrp);
287 gmx_fio_do_int(fio,rot->nstrout);
288 gmx_fio_do_int(fio,rot->nstsout);
290 snew(rot->grp,rot->ngrp);
291 for(g=0; g<rot->ngrp; g++)
292 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
296 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
297 int file_version, real *fudgeQQ)
299 int i,j,k,*tmp,idum=0;
304 real zerotemptime,finish_t,init_temp,finish_temp;
306 if (file_version != tpx_version)
308 /* Give a warning about features that are not accessible */
309 fprintf(stderr,"Note: tpx file_version %d, software version %d\n",
310 file_version,tpx_version);
318 if (file_version == 0)
323 /* Basic inputrec stuff */
324 gmx_fio_do_int(fio,ir->eI);
325 if (file_version >= 62) {
326 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
328 gmx_fio_do_int(fio,idum);
331 if(file_version > 25) {
332 if (file_version >= 62) {
333 gmx_fio_do_gmx_large_int(fio, ir->init_step);
335 gmx_fio_do_int(fio,idum);
336 ir->init_step = idum;
342 if(file_version >= 58)
343 gmx_fio_do_int(fio,ir->simulation_part);
345 ir->simulation_part=1;
347 if (file_version >= 67) {
348 gmx_fio_do_int(fio,ir->nstcalcenergy);
350 ir->nstcalcenergy = 1;
352 if (file_version < 53) {
353 /* The pbc info has been moved out of do_inputrec,
354 * since we always want it, also without reading the inputrec.
356 gmx_fio_do_int(fio,ir->ePBC);
357 if ((file_version <= 15) && (ir->ePBC == 2))
359 if (file_version >= 45) {
360 gmx_fio_do_int(fio,ir->bPeriodicMols);
364 ir->bPeriodicMols = TRUE;
366 ir->bPeriodicMols = FALSE;
370 gmx_fio_do_int(fio,ir->ns_type);
371 gmx_fio_do_int(fio,ir->nstlist);
372 gmx_fio_do_int(fio,ir->ndelta);
373 if (file_version < 41) {
374 gmx_fio_do_int(fio,idum);
375 gmx_fio_do_int(fio,idum);
377 if (file_version >= 45)
378 gmx_fio_do_real(fio,ir->rtpi);
381 gmx_fio_do_int(fio,ir->nstcomm);
382 if (file_version > 34)
383 gmx_fio_do_int(fio,ir->comm_mode);
384 else if (ir->nstcomm < 0)
385 ir->comm_mode = ecmANGULAR;
387 ir->comm_mode = ecmLINEAR;
388 ir->nstcomm = abs(ir->nstcomm);
390 if(file_version > 25)
391 gmx_fio_do_int(fio,ir->nstcheckpoint);
395 gmx_fio_do_int(fio,ir->nstcgsteep);
398 gmx_fio_do_int(fio,ir->nbfgscorr);
402 gmx_fio_do_int(fio,ir->nstlog);
403 gmx_fio_do_int(fio,ir->nstxout);
404 gmx_fio_do_int(fio,ir->nstvout);
405 gmx_fio_do_int(fio,ir->nstfout);
406 gmx_fio_do_int(fio,ir->nstenergy);
407 gmx_fio_do_int(fio,ir->nstxtcout);
408 if (file_version >= 59) {
409 gmx_fio_do_double(fio,ir->init_t);
410 gmx_fio_do_double(fio,ir->delta_t);
412 gmx_fio_do_real(fio,rdum);
414 gmx_fio_do_real(fio,rdum);
417 gmx_fio_do_real(fio,ir->xtcprec);
418 if (file_version < 19) {
419 gmx_fio_do_int(fio,idum);
420 gmx_fio_do_int(fio,idum);
422 if(file_version < 18)
423 gmx_fio_do_int(fio,idum);
424 gmx_fio_do_real(fio,ir->rlist);
425 if (file_version >= 67) {
426 gmx_fio_do_real(fio,ir->rlistlong);
428 gmx_fio_do_int(fio,ir->coulombtype);
429 if (file_version < 32 && ir->coulombtype == eelRF)
430 ir->coulombtype = eelRF_NEC;
431 gmx_fio_do_real(fio,ir->rcoulomb_switch);
432 gmx_fio_do_real(fio,ir->rcoulomb);
433 gmx_fio_do_int(fio,ir->vdwtype);
434 gmx_fio_do_real(fio,ir->rvdw_switch);
435 gmx_fio_do_real(fio,ir->rvdw);
436 if (file_version < 67) {
437 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
439 gmx_fio_do_int(fio,ir->eDispCorr);
440 gmx_fio_do_real(fio,ir->epsilon_r);
441 if (file_version >= 37) {
442 gmx_fio_do_real(fio,ir->epsilon_rf);
444 if (EEL_RF(ir->coulombtype)) {
445 ir->epsilon_rf = ir->epsilon_r;
448 ir->epsilon_rf = 1.0;
451 if (file_version >= 29)
452 gmx_fio_do_real(fio,ir->tabext);
456 if(file_version > 25) {
457 gmx_fio_do_int(fio,ir->gb_algorithm);
458 gmx_fio_do_int(fio,ir->nstgbradii);
459 gmx_fio_do_real(fio,ir->rgbradii);
460 gmx_fio_do_real(fio,ir->gb_saltconc);
461 gmx_fio_do_int(fio,ir->implicit_solvent);
463 ir->gb_algorithm=egbSTILL;
467 ir->implicit_solvent=eisNO;
471 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
472 gmx_fio_do_real(fio,ir->gb_obc_alpha);
473 gmx_fio_do_real(fio,ir->gb_obc_beta);
474 gmx_fio_do_real(fio,ir->gb_obc_gamma);
477 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
478 gmx_fio_do_int(fio,ir->sa_algorithm);
482 ir->gb_dielectric_offset = 0.009;
483 ir->sa_algorithm = esaAPPROX;
485 gmx_fio_do_real(fio,ir->sa_surface_tension);
487 /* Override sa_surface_tension if it is not changed in the mpd-file */
488 if(ir->sa_surface_tension<0)
490 if(ir->gb_algorithm==egbSTILL)
492 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
494 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
496 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
503 /* Better use sensible values than insane (0.0) ones... */
504 ir->gb_epsilon_solvent = 80;
505 ir->gb_obc_alpha = 1.0;
506 ir->gb_obc_beta = 0.8;
507 ir->gb_obc_gamma = 4.85;
508 ir->sa_surface_tension = 2.092;
512 gmx_fio_do_int(fio,ir->nkx);
513 gmx_fio_do_int(fio,ir->nky);
514 gmx_fio_do_int(fio,ir->nkz);
515 gmx_fio_do_int(fio,ir->pme_order);
516 gmx_fio_do_real(fio,ir->ewald_rtol);
518 if (file_version >=24)
519 gmx_fio_do_int(fio,ir->ewald_geometry);
521 ir->ewald_geometry=eewg3D;
523 if (file_version <=17) {
524 ir->epsilon_surface=0;
525 if (file_version==17)
526 gmx_fio_do_int(fio,idum);
529 gmx_fio_do_real(fio,ir->epsilon_surface);
531 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
533 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
534 gmx_fio_do_int(fio,ir->etc);
535 /* before version 18, ir->etc was a gmx_bool (ir->btc),
536 * but the values 0 and 1 still mean no and
537 * berendsen temperature coupling, respectively.
539 if (file_version >= 71)
541 gmx_fio_do_int(fio,ir->nsttcouple);
545 ir->nsttcouple = ir->nstcalcenergy;
547 if (file_version <= 15)
549 gmx_fio_do_int(fio,idum);
551 if (file_version <=17)
553 gmx_fio_do_int(fio,ir->epct);
554 if (file_version <= 15)
558 ir->epct = epctSURFACETENSION;
560 gmx_fio_do_int(fio,idum);
563 /* we have removed the NO alternative at the beginning */
567 ir->epct=epctISOTROPIC;
571 ir->epc=epcBERENDSEN;
576 gmx_fio_do_int(fio,ir->epc);
577 gmx_fio_do_int(fio,ir->epct);
579 if (file_version >= 71)
581 gmx_fio_do_int(fio,ir->nstpcouple);
585 ir->nstpcouple = ir->nstcalcenergy;
587 gmx_fio_do_real(fio,ir->tau_p);
588 if (file_version <= 15) {
589 gmx_fio_do_rvec(fio,vdum);
590 clear_mat(ir->ref_p);
592 ir->ref_p[i][i] = vdum[i];
594 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
595 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
596 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
598 if (file_version <= 15) {
599 gmx_fio_do_rvec(fio,vdum);
600 clear_mat(ir->compress);
602 ir->compress[i][i] = vdum[i];
605 gmx_fio_do_rvec(fio,ir->compress[XX]);
606 gmx_fio_do_rvec(fio,ir->compress[YY]);
607 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
609 if (file_version >= 47) {
610 gmx_fio_do_int(fio,ir->refcoord_scaling);
611 gmx_fio_do_rvec(fio,ir->posres_com);
612 gmx_fio_do_rvec(fio,ir->posres_comB);
614 ir->refcoord_scaling = erscNO;
615 clear_rvec(ir->posres_com);
616 clear_rvec(ir->posres_comB);
618 if(file_version > 25)
619 gmx_fio_do_int(fio,ir->andersen_seed);
623 if(file_version < 26) {
624 gmx_fio_do_gmx_bool(fio,bSimAnn);
625 gmx_fio_do_real(fio,zerotemptime);
628 if (file_version < 37)
629 gmx_fio_do_real(fio,rdum);
631 gmx_fio_do_real(fio,ir->shake_tol);
632 if (file_version < 54)
633 gmx_fio_do_real(fio,*fudgeQQ);
634 gmx_fio_do_int(fio,ir->efep);
635 if (file_version <= 14 && ir->efep > efepNO)
637 if (file_version >= 59) {
638 gmx_fio_do_double(fio, ir->init_lambda);
639 gmx_fio_do_double(fio, ir->delta_lambda);
641 gmx_fio_do_real(fio,rdum);
642 ir->init_lambda = rdum;
643 gmx_fio_do_real(fio,rdum);
644 ir->delta_lambda = rdum;
646 if (file_version >= 64) {
647 gmx_fio_do_int(fio,ir->n_flambda);
649 snew(ir->flambda,ir->n_flambda);
651 bDum=gmx_fio_ndo_double(fio,ir->flambda,ir->n_flambda);
656 if (file_version >= 13)
657 gmx_fio_do_real(fio,ir->sc_alpha);
660 if (file_version >= 38)
661 gmx_fio_do_int(fio,ir->sc_power);
664 if (file_version >= 15)
665 gmx_fio_do_real(fio,ir->sc_sigma);
670 if (file_version >= 71)
672 ir->sc_sigma_min = ir->sc_sigma;
676 ir->sc_sigma_min = 0;
679 if (file_version >= 64) {
680 gmx_fio_do_int(fio,ir->nstdhdl);
685 if (file_version >= 73)
687 gmx_fio_do_int(fio, ir->separate_dhdl_file);
688 gmx_fio_do_int(fio, ir->dhdl_derivatives);
692 ir->separate_dhdl_file = sepdhdlfileYES;
693 ir->dhdl_derivatives = dhdlderivativesYES;
696 if (file_version >= 71)
698 gmx_fio_do_int(fio,ir->dh_hist_size);
699 gmx_fio_do_double(fio,ir->dh_hist_spacing);
703 ir->dh_hist_size = 0;
704 ir->dh_hist_spacing = 0.1;
706 if (file_version >= 57) {
707 gmx_fio_do_int(fio,ir->eDisre);
709 gmx_fio_do_int(fio,ir->eDisreWeighting);
710 if (file_version < 22) {
711 if (ir->eDisreWeighting == 0)
712 ir->eDisreWeighting = edrwEqual;
714 ir->eDisreWeighting = edrwConservative;
716 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
717 gmx_fio_do_real(fio,ir->dr_fc);
718 gmx_fio_do_real(fio,ir->dr_tau);
719 gmx_fio_do_int(fio,ir->nstdisreout);
720 if (file_version >= 22) {
721 gmx_fio_do_real(fio,ir->orires_fc);
722 gmx_fio_do_real(fio,ir->orires_tau);
723 gmx_fio_do_int(fio,ir->nstorireout);
729 if(file_version >= 26) {
730 gmx_fio_do_real(fio,ir->dihre_fc);
731 if (file_version < 56) {
732 gmx_fio_do_real(fio,rdum);
733 gmx_fio_do_int(fio,idum);
739 gmx_fio_do_real(fio,ir->em_stepsize);
740 gmx_fio_do_real(fio,ir->em_tol);
741 if (file_version >= 22)
742 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
744 ir->bShakeSOR = TRUE;
745 if (file_version >= 11)
746 gmx_fio_do_int(fio,ir->niter);
749 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
752 if (file_version >= 21)
753 gmx_fio_do_real(fio,ir->fc_stepsize);
756 gmx_fio_do_int(fio,ir->eConstrAlg);
757 gmx_fio_do_int(fio,ir->nProjOrder);
758 gmx_fio_do_real(fio,ir->LincsWarnAngle);
759 if (file_version <= 14)
760 gmx_fio_do_int(fio,idum);
761 if (file_version >=26)
762 gmx_fio_do_int(fio,ir->nLincsIter);
765 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
768 if (file_version < 33)
769 gmx_fio_do_real(fio,bd_temp);
770 gmx_fio_do_real(fio,ir->bd_fric);
771 gmx_fio_do_int(fio,ir->ld_seed);
772 if (file_version >= 33) {
774 gmx_fio_do_rvec(fio,ir->deform[i]);
777 clear_rvec(ir->deform[i]);
779 if (file_version >= 14)
780 gmx_fio_do_real(fio,ir->cos_accel);
783 gmx_fio_do_int(fio,ir->userint1);
784 gmx_fio_do_int(fio,ir->userint2);
785 gmx_fio_do_int(fio,ir->userint3);
786 gmx_fio_do_int(fio,ir->userint4);
787 gmx_fio_do_real(fio,ir->userreal1);
788 gmx_fio_do_real(fio,ir->userreal2);
789 gmx_fio_do_real(fio,ir->userreal3);
790 gmx_fio_do_real(fio,ir->userreal4);
793 if (file_version >= 48) {
794 gmx_fio_do_int(fio,ir->ePull);
795 if (ir->ePull != epullNO) {
798 do_pull(fio, ir->pull,bRead,file_version);
804 /* Enforced rotation */
805 if (file_version >= 74) {
806 gmx_fio_do_int(fio,ir->bRot);
807 if (ir->bRot == TRUE) {
810 do_rot(fio, ir->rot,bRead,file_version);
817 gmx_fio_do_int(fio,ir->opts.ngtc);
818 if (file_version >= 69) {
819 gmx_fio_do_int(fio,ir->opts.nhchainlength);
821 ir->opts.nhchainlength = 1;
823 gmx_fio_do_int(fio,ir->opts.ngacc);
824 gmx_fio_do_int(fio,ir->opts.ngfrz);
825 gmx_fio_do_int(fio,ir->opts.ngener);
828 snew(ir->opts.nrdf, ir->opts.ngtc);
829 snew(ir->opts.ref_t, ir->opts.ngtc);
830 snew(ir->opts.annealing, ir->opts.ngtc);
831 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
832 snew(ir->opts.anneal_time, ir->opts.ngtc);
833 snew(ir->opts.anneal_temp, ir->opts.ngtc);
834 snew(ir->opts.tau_t, ir->opts.ngtc);
835 snew(ir->opts.nFreeze,ir->opts.ngfrz);
836 snew(ir->opts.acc, ir->opts.ngacc);
837 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
839 if (ir->opts.ngtc > 0) {
840 if (bRead && file_version<13) {
841 snew(tmp,ir->opts.ngtc);
842 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
843 for(i=0; i<ir->opts.ngtc; i++)
844 ir->opts.nrdf[i] = tmp[i];
847 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
849 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
850 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
851 if (file_version<33 && ir->eI==eiBD) {
852 for(i=0; i<ir->opts.ngtc; i++)
853 ir->opts.tau_t[i] = bd_temp;
856 if (ir->opts.ngfrz > 0)
857 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
858 if (ir->opts.ngacc > 0)
859 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
860 if (file_version >= 12)
861 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
862 ir->opts.ngener*ir->opts.ngener);
864 if(bRead && file_version < 26) {
865 for(i=0;i<ir->opts.ngtc;i++) {
867 ir->opts.annealing[i] = eannSINGLE;
868 ir->opts.anneal_npoints[i] = 2;
869 snew(ir->opts.anneal_time[i],2);
870 snew(ir->opts.anneal_temp[i],2);
871 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
872 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
873 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
874 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
875 ir->opts.anneal_time[i][0] = ir->init_t;
876 ir->opts.anneal_time[i][1] = finish_t;
877 ir->opts.anneal_temp[i][0] = init_temp;
878 ir->opts.anneal_temp[i][1] = finish_temp;
880 ir->opts.annealing[i] = eannNO;
881 ir->opts.anneal_npoints[i] = 0;
885 /* file version 26 or later */
886 /* First read the lists with annealing and npoints for each group */
887 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
888 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
889 for(j=0;j<(ir->opts.ngtc);j++) {
890 k=ir->opts.anneal_npoints[j];
892 snew(ir->opts.anneal_time[j],k);
893 snew(ir->opts.anneal_temp[j],k);
895 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
896 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
900 if (file_version >= 45) {
901 gmx_fio_do_int(fio,ir->nwall);
902 gmx_fio_do_int(fio,ir->wall_type);
903 if (file_version >= 50)
904 gmx_fio_do_real(fio,ir->wall_r_linpot);
906 ir->wall_r_linpot = -1;
907 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
908 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
909 gmx_fio_do_real(fio,ir->wall_density[0]);
910 gmx_fio_do_real(fio,ir->wall_density[1]);
911 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
915 ir->wall_atomtype[0] = -1;
916 ir->wall_atomtype[1] = -1;
917 ir->wall_density[0] = 0;
918 ir->wall_density[1] = 0;
919 ir->wall_ewald_zfac = 3;
921 /* Cosine stuff for electric fields */
922 for(j=0; (j<DIM); j++) {
923 gmx_fio_do_int(fio,ir->ex[j].n);
924 gmx_fio_do_int(fio,ir->et[j].n);
926 snew(ir->ex[j].a, ir->ex[j].n);
927 snew(ir->ex[j].phi,ir->ex[j].n);
928 snew(ir->et[j].a, ir->et[j].n);
929 snew(ir->et[j].phi,ir->et[j].n);
931 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
932 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
933 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
934 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
938 if(file_version>=39){
939 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
940 gmx_fio_do_int(fio,ir->QMMMscheme);
941 gmx_fio_do_real(fio,ir->scalefactor);
942 gmx_fio_do_int(fio,ir->opts.ngQM);
944 snew(ir->opts.QMmethod, ir->opts.ngQM);
945 snew(ir->opts.QMbasis, ir->opts.ngQM);
946 snew(ir->opts.QMcharge, ir->opts.ngQM);
947 snew(ir->opts.QMmult, ir->opts.ngQM);
948 snew(ir->opts.bSH, ir->opts.ngQM);
949 snew(ir->opts.CASorbitals, ir->opts.ngQM);
950 snew(ir->opts.CASelectrons,ir->opts.ngQM);
951 snew(ir->opts.SAon, ir->opts.ngQM);
952 snew(ir->opts.SAoff, ir->opts.ngQM);
953 snew(ir->opts.SAsteps, ir->opts.ngQM);
954 snew(ir->opts.bOPT, ir->opts.ngQM);
955 snew(ir->opts.bTS, ir->opts.ngQM);
957 if (ir->opts.ngQM > 0) {
958 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
959 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
960 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
961 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
962 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
963 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
964 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
965 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
966 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
967 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
968 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
969 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
971 /* end of QMMM stuff */
976 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
978 gmx_fio_do_real(fio,iparams->harmonic.rA);
979 gmx_fio_do_real(fio,iparams->harmonic.krA);
980 gmx_fio_do_real(fio,iparams->harmonic.rB);
981 gmx_fio_do_real(fio,iparams->harmonic.krB);
984 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
985 gmx_bool bRead, int file_version)
992 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1000 do_harm(fio, iparams,bRead);
1001 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1002 /* Correct incorrect storage of parameters */
1003 iparams->pdihs.phiB = iparams->pdihs.phiA;
1004 iparams->pdihs.cpB = iparams->pdihs.cpA;
1008 gmx_fio_do_real(fio,iparams->fene.bm);
1009 gmx_fio_do_real(fio,iparams->fene.kb);
1012 gmx_fio_do_real(fio,iparams->restraint.lowA);
1013 gmx_fio_do_real(fio,iparams->restraint.up1A);
1014 gmx_fio_do_real(fio,iparams->restraint.up2A);
1015 gmx_fio_do_real(fio,iparams->restraint.kA);
1016 gmx_fio_do_real(fio,iparams->restraint.lowB);
1017 gmx_fio_do_real(fio,iparams->restraint.up1B);
1018 gmx_fio_do_real(fio,iparams->restraint.up2B);
1019 gmx_fio_do_real(fio,iparams->restraint.kB);
1025 gmx_fio_do_real(fio,iparams->tab.kA);
1026 gmx_fio_do_int(fio,iparams->tab.table);
1027 gmx_fio_do_real(fio,iparams->tab.kB);
1029 case F_CROSS_BOND_BONDS:
1030 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1031 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1032 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1034 case F_CROSS_BOND_ANGLES:
1035 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1036 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1037 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1038 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1040 case F_UREY_BRADLEY:
1041 gmx_fio_do_real(fio,iparams->u_b.theta);
1042 gmx_fio_do_real(fio,iparams->u_b.ktheta);
1043 gmx_fio_do_real(fio,iparams->u_b.r13);
1044 gmx_fio_do_real(fio,iparams->u_b.kUB);
1046 case F_QUARTIC_ANGLES:
1047 gmx_fio_do_real(fio,iparams->qangle.theta);
1048 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1051 gmx_fio_do_real(fio,iparams->bham.a);
1052 gmx_fio_do_real(fio,iparams->bham.b);
1053 gmx_fio_do_real(fio,iparams->bham.c);
1056 gmx_fio_do_real(fio,iparams->morse.b0);
1057 gmx_fio_do_real(fio,iparams->morse.cb);
1058 gmx_fio_do_real(fio,iparams->morse.beta);
1061 gmx_fio_do_real(fio,iparams->cubic.b0);
1062 gmx_fio_do_real(fio,iparams->cubic.kb);
1063 gmx_fio_do_real(fio,iparams->cubic.kcub);
1067 case F_POLARIZATION:
1068 gmx_fio_do_real(fio,iparams->polarize.alpha);
1071 if (file_version < 31)
1072 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1073 gmx_fio_do_real(fio,iparams->wpol.al_x);
1074 gmx_fio_do_real(fio,iparams->wpol.al_y);
1075 gmx_fio_do_real(fio,iparams->wpol.al_z);
1076 gmx_fio_do_real(fio,iparams->wpol.rOH);
1077 gmx_fio_do_real(fio,iparams->wpol.rHH);
1078 gmx_fio_do_real(fio,iparams->wpol.rOD);
1081 gmx_fio_do_real(fio,iparams->thole.a);
1082 gmx_fio_do_real(fio,iparams->thole.alpha1);
1083 gmx_fio_do_real(fio,iparams->thole.alpha2);
1084 gmx_fio_do_real(fio,iparams->thole.rfac);
1087 gmx_fio_do_real(fio,iparams->lj.c6);
1088 gmx_fio_do_real(fio,iparams->lj.c12);
1091 gmx_fio_do_real(fio,iparams->lj14.c6A);
1092 gmx_fio_do_real(fio,iparams->lj14.c12A);
1093 gmx_fio_do_real(fio,iparams->lj14.c6B);
1094 gmx_fio_do_real(fio,iparams->lj14.c12B);
1097 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1098 gmx_fio_do_real(fio,iparams->ljc14.qi);
1099 gmx_fio_do_real(fio,iparams->ljc14.qj);
1100 gmx_fio_do_real(fio,iparams->ljc14.c6);
1101 gmx_fio_do_real(fio,iparams->ljc14.c12);
1103 case F_LJC_PAIRS_NB:
1104 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1105 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1106 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1107 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1113 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1114 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1115 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1116 /* Read the incorrectly stored multiplicity */
1117 gmx_fio_do_real(fio,iparams->harmonic.rB);
1118 gmx_fio_do_real(fio,iparams->harmonic.krB);
1119 iparams->pdihs.phiB = iparams->pdihs.phiA;
1120 iparams->pdihs.cpB = iparams->pdihs.cpA;
1122 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1123 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1124 gmx_fio_do_int(fio,iparams->pdihs.mult);
1128 gmx_fio_do_int(fio,iparams->disres.label);
1129 gmx_fio_do_int(fio,iparams->disres.type);
1130 gmx_fio_do_real(fio,iparams->disres.low);
1131 gmx_fio_do_real(fio,iparams->disres.up1);
1132 gmx_fio_do_real(fio,iparams->disres.up2);
1133 gmx_fio_do_real(fio,iparams->disres.kfac);
1136 gmx_fio_do_int(fio,iparams->orires.ex);
1137 gmx_fio_do_int(fio,iparams->orires.label);
1138 gmx_fio_do_int(fio,iparams->orires.power);
1139 gmx_fio_do_real(fio,iparams->orires.c);
1140 gmx_fio_do_real(fio,iparams->orires.obs);
1141 gmx_fio_do_real(fio,iparams->orires.kfac);
1144 gmx_fio_do_int(fio,iparams->dihres.power);
1145 gmx_fio_do_int(fio,iparams->dihres.label);
1146 gmx_fio_do_real(fio,iparams->dihres.phi);
1147 gmx_fio_do_real(fio,iparams->dihres.dphi);
1148 gmx_fio_do_real(fio,iparams->dihres.kfac);
1151 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1152 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1153 if (bRead && file_version < 27) {
1154 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1155 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1157 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1158 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1162 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1163 if(file_version>=25)
1164 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1167 /* Fourier dihedrals are internally represented
1168 * as Ryckaert-Bellemans since those are faster to compute.
1170 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1171 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1175 gmx_fio_do_real(fio,iparams->constr.dA);
1176 gmx_fio_do_real(fio,iparams->constr.dB);
1179 gmx_fio_do_real(fio,iparams->settle.doh);
1180 gmx_fio_do_real(fio,iparams->settle.dhh);
1183 gmx_fio_do_real(fio,iparams->vsite.a);
1188 gmx_fio_do_real(fio,iparams->vsite.a);
1189 gmx_fio_do_real(fio,iparams->vsite.b);
1194 gmx_fio_do_real(fio,iparams->vsite.a);
1195 gmx_fio_do_real(fio,iparams->vsite.b);
1196 gmx_fio_do_real(fio,iparams->vsite.c);
1199 gmx_fio_do_int(fio,iparams->vsiten.n);
1200 gmx_fio_do_real(fio,iparams->vsiten.a);
1205 /* We got rid of some parameters in version 68 */
1206 if(bRead && file_version<68)
1208 gmx_fio_do_real(fio,rdum);
1209 gmx_fio_do_real(fio,rdum);
1210 gmx_fio_do_real(fio,rdum);
1211 gmx_fio_do_real(fio,rdum);
1213 gmx_fio_do_real(fio,iparams->gb.sar);
1214 gmx_fio_do_real(fio,iparams->gb.st);
1215 gmx_fio_do_real(fio,iparams->gb.pi);
1216 gmx_fio_do_real(fio,iparams->gb.gbr);
1217 gmx_fio_do_real(fio,iparams->gb.bmlt);
1220 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1221 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1224 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1226 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1229 gmx_fio_unset_comment(fio);
1232 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1239 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1241 if (file_version < 44) {
1242 for(i=0; i<MAXNODES; i++)
1243 gmx_fio_do_int(fio,idum);
1245 gmx_fio_do_int(fio,ilist->nr);
1247 snew(ilist->iatoms,ilist->nr);
1248 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1250 gmx_fio_unset_comment(fio);
1253 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1254 gmx_bool bRead, int file_version)
1260 gmx_fio_do_int(fio,ffparams->atnr);
1261 if (file_version < 57) {
1262 gmx_fio_do_int(fio,idum);
1264 gmx_fio_do_int(fio,ffparams->ntypes);
1266 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1267 ffparams->atnr,ffparams->ntypes);
1269 snew(ffparams->functype,ffparams->ntypes);
1270 snew(ffparams->iparams,ffparams->ntypes);
1272 /* Read/write all the function types */
1273 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1275 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1277 if (file_version >= 66) {
1278 gmx_fio_do_double(fio,ffparams->reppow);
1280 ffparams->reppow = 12.0;
1283 if (file_version >= 57) {
1284 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1287 /* Check whether all these function types are supported by the code.
1288 * In practice the code is backwards compatible, which means that the
1289 * numbering may have to be altered from old numbering to new numbering
1291 for (i=0; (i<ffparams->ntypes); i++) {
1293 /* Loop over file versions */
1294 for (k=0; (k<NFTUPD); k++)
1295 /* Compare the read file_version to the update table */
1296 if ((file_version < ftupd[k].fvnr) &&
1297 (ffparams->functype[i] >= ftupd[k].ftype)) {
1298 ffparams->functype[i] += 1;
1300 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1301 i,ffparams->functype[i],
1302 interaction_function[ftupd[k].ftype].longname);
1307 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1310 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1314 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1317 int i,j,renum[F_NRE];
1318 gmx_bool bDum=TRUE,bClear;
1321 for(j=0; (j<F_NRE); j++) {
1324 for (k=0; k<NFTUPD; k++)
1325 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1329 ilist[j].iatoms = NULL;
1331 do_ilist(fio, &ilist[j],bRead,file_version,j);
1334 if (bRead && gmx_debug_at)
1335 pr_ilist(debug,0,interaction_function[j].longname,
1336 functype,&ilist[j],TRUE);
1341 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1342 gmx_bool bRead, int file_version)
1344 do_ffparams(fio, ffparams,bRead,file_version);
1346 if (file_version >= 54) {
1347 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1350 do_ilists(fio, molt->ilist,bRead,file_version);
1353 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1355 int i,idum,dum_nra,*dum_a;
1358 if (file_version < 44)
1359 for(i=0; i<MAXNODES; i++)
1360 gmx_fio_do_int(fio,idum);
1361 gmx_fio_do_int(fio,block->nr);
1362 if (file_version < 51)
1363 gmx_fio_do_int(fio,dum_nra);
1365 block->nalloc_index = block->nr+1;
1366 snew(block->index,block->nalloc_index);
1368 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1370 if (file_version < 51 && dum_nra > 0) {
1371 snew(dum_a,dum_nra);
1372 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1377 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1383 if (file_version < 44)
1384 for(i=0; i<MAXNODES; i++)
1385 gmx_fio_do_int(fio,idum);
1386 gmx_fio_do_int(fio,block->nr);
1387 gmx_fio_do_int(fio,block->nra);
1389 block->nalloc_index = block->nr+1;
1390 snew(block->index,block->nalloc_index);
1391 block->nalloc_a = block->nra;
1392 snew(block->a,block->nalloc_a);
1394 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1395 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1398 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1399 int file_version, gmx_groups_t *groups,int atnr)
1403 gmx_fio_do_real(fio,atom->m);
1404 gmx_fio_do_real(fio,atom->q);
1405 gmx_fio_do_real(fio,atom->mB);
1406 gmx_fio_do_real(fio,atom->qB);
1407 gmx_fio_do_ushort(fio, atom->type);
1408 gmx_fio_do_ushort(fio, atom->typeB);
1409 gmx_fio_do_int(fio,atom->ptype);
1410 gmx_fio_do_int(fio,atom->resind);
1411 if (file_version >= 52)
1412 gmx_fio_do_int(fio,atom->atomnumber);
1414 atom->atomnumber = NOTSET;
1415 if (file_version < 23)
1417 else if (file_version < 39)
1422 if (file_version < 57) {
1423 unsigned char uchar[egcNR];
1424 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1425 for(i=myngrp; (i<ngrp); i++) {
1428 /* Copy the old data format to the groups struct */
1429 for(i=0; i<ngrp; i++) {
1430 groups->grpnr[i][atnr] = uchar[i];
1435 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1441 if (file_version < 23)
1443 else if (file_version < 39)
1448 for(j=0; (j<ngrp); j++) {
1450 gmx_fio_do_int(fio,grps[j].nr);
1452 snew(grps[j].nm_ind,grps[j].nr);
1453 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1457 snew(grps[j].nm_ind,grps[j].nr);
1462 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1467 gmx_fio_do_int(fio,ls);
1468 *nm = get_symtab_handle(symtab,ls);
1471 ls = lookup_symtab(symtab,*nm);
1472 gmx_fio_do_int(fio,ls);
1476 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1481 for (j=0; (j<nstr); j++)
1482 do_symstr(fio, &(nm[j]),bRead,symtab);
1485 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1486 t_symtab *symtab, int file_version)
1490 for (j=0; (j<n); j++) {
1491 do_symstr(fio, &(ri[j].name),bRead,symtab);
1492 if (file_version >= 63) {
1493 gmx_fio_do_int(fio,ri[j].nr);
1494 gmx_fio_do_uchar(fio, ri[j].ic);
1502 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1504 gmx_groups_t *groups)
1508 gmx_fio_do_int(fio,atoms->nr);
1509 gmx_fio_do_int(fio,atoms->nres);
1510 if (file_version < 57) {
1511 gmx_fio_do_int(fio,groups->ngrpname);
1512 for(i=0; i<egcNR; i++) {
1513 groups->ngrpnr[i] = atoms->nr;
1514 snew(groups->grpnr[i],groups->ngrpnr[i]);
1518 snew(atoms->atom,atoms->nr);
1519 snew(atoms->atomname,atoms->nr);
1520 snew(atoms->atomtype,atoms->nr);
1521 snew(atoms->atomtypeB,atoms->nr);
1522 snew(atoms->resinfo,atoms->nres);
1523 if (file_version < 57) {
1524 snew(groups->grpname,groups->ngrpname);
1526 atoms->pdbinfo = NULL;
1528 for(i=0; (i<atoms->nr); i++) {
1529 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1531 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1532 if (bRead && (file_version <= 20)) {
1533 for(i=0; i<atoms->nr; i++) {
1534 atoms->atomtype[i] = put_symtab(symtab,"?");
1535 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1538 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1539 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1541 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1543 if (file_version < 57) {
1544 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1546 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1550 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1551 gmx_bool bRead,t_symtab *symtab,
1557 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1558 gmx_fio_do_int(fio,groups->ngrpname);
1560 snew(groups->grpname,groups->ngrpname);
1562 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1563 for(g=0; g<egcNR; g++) {
1564 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1565 if (groups->ngrpnr[g] == 0) {
1567 groups->grpnr[g] = NULL;
1571 snew(groups->grpnr[g],groups->ngrpnr[g]);
1573 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1578 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1579 t_symtab *symtab,int file_version)
1582 gmx_bool bDum = TRUE;
1584 if (file_version > 25) {
1585 gmx_fio_do_int(fio,atomtypes->nr);
1588 snew(atomtypes->radius,j);
1589 snew(atomtypes->vol,j);
1590 snew(atomtypes->surftens,j);
1591 snew(atomtypes->atomnumber,j);
1592 snew(atomtypes->gb_radius,j);
1593 snew(atomtypes->S_hct,j);
1595 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1596 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1597 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1598 if(file_version >= 40)
1600 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1602 if(file_version >= 60)
1604 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1605 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1608 /* File versions prior to 26 cannot do GBSA,
1609 * so they dont use this structure
1612 atomtypes->radius = NULL;
1613 atomtypes->vol = NULL;
1614 atomtypes->surftens = NULL;
1615 atomtypes->atomnumber = NULL;
1616 atomtypes->gb_radius = NULL;
1617 atomtypes->S_hct = NULL;
1621 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1627 gmx_fio_do_int(fio,symtab->nr);
1630 snew(symtab->symbuf,1);
1631 symbuf = symtab->symbuf;
1632 symbuf->bufsize = nr;
1633 snew(symbuf->buf,nr);
1634 for (i=0; (i<nr); i++) {
1635 gmx_fio_do_string(fio,buf);
1636 symbuf->buf[i]=strdup(buf);
1640 symbuf = symtab->symbuf;
1641 while (symbuf!=NULL) {
1642 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1643 gmx_fio_do_string(fio,symbuf->buf[i]);
1645 symbuf=symbuf->next;
1648 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1652 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1654 int i,j,ngrid,gs,nelem;
1656 gmx_fio_do_int(fio,cmap_grid->ngrid);
1657 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1659 ngrid = cmap_grid->ngrid;
1660 gs = cmap_grid->grid_spacing;
1665 snew(cmap_grid->cmapdata,ngrid);
1667 for(i=0;i<cmap_grid->ngrid;i++)
1669 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1673 for(i=0;i<cmap_grid->ngrid;i++)
1675 for(j=0;j<nelem;j++)
1677 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1678 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1679 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1680 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1686 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1692 /* We always assign a new chain number, but save the chain id characters
1693 * for larger molecules.
1695 #define CHAIN_MIN_ATOMS 15
1699 for(m=0; m<mols->nr; m++)
1702 a1=mols->index[m+1];
1703 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
1712 for(a=a0; a<a1; a++)
1714 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
1715 atoms->resinfo[atoms->atom[a].resind].chainid = c;
1720 /* Blank out the chain id if there was only one chain */
1723 for(r=0; r<atoms->nres; r++)
1725 atoms->resinfo[r].chainid = ' ';
1730 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
1731 t_symtab *symtab, int file_version,
1732 gmx_groups_t *groups)
1736 if (file_version >= 57) {
1737 do_symstr(fio, &(molt->name),bRead,symtab);
1740 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
1742 if (bRead && gmx_debug_at) {
1743 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
1746 if (file_version >= 57) {
1747 do_ilists(fio, molt->ilist,bRead,file_version);
1749 do_block(fio, &molt->cgs,bRead,file_version);
1750 if (bRead && gmx_debug_at) {
1751 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
1755 /* This used to be in the atoms struct */
1756 do_blocka(fio, &molt->excls, bRead, file_version);
1759 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
1764 gmx_fio_do_int(fio,molb->type);
1765 gmx_fio_do_int(fio,molb->nmol);
1766 gmx_fio_do_int(fio,molb->natoms_mol);
1767 /* Position restraint coordinates */
1768 gmx_fio_do_int(fio,molb->nposres_xA);
1769 if (molb->nposres_xA > 0) {
1771 snew(molb->posres_xA,molb->nposres_xA);
1773 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
1775 gmx_fio_do_int(fio,molb->nposres_xB);
1776 if (molb->nposres_xB > 0) {
1778 snew(molb->posres_xB,molb->nposres_xB);
1780 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
1785 static t_block mtop_mols(gmx_mtop_t *mtop)
1791 for(mb=0; mb<mtop->nmolblock; mb++) {
1792 mols.nr += mtop->molblock[mb].nmol;
1794 mols.nalloc_index = mols.nr + 1;
1795 snew(mols.index,mols.nalloc_index);
1800 for(mb=0; mb<mtop->nmolblock; mb++) {
1801 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
1802 a += mtop->molblock[mb].natoms_mol;
1811 static void add_posres_molblock(gmx_mtop_t *mtop)
1816 gmx_molblock_t *molb;
1819 il = &mtop->moltype[0].ilist[F_POSRES];
1825 for(i=0; i<il->nr; i+=2) {
1826 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1827 am = max(am,il->iatoms[i+1]);
1828 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
1829 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
1830 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
1834 /* Make the posres coordinate block end at a molecule end */
1836 while(am >= mtop->mols.index[mol+1]) {
1839 molb = &mtop->molblock[0];
1840 molb->nposres_xA = mtop->mols.index[mol+1];
1841 snew(molb->posres_xA,molb->nposres_xA);
1843 molb->nposres_xB = molb->nposres_xA;
1844 snew(molb->posres_xB,molb->nposres_xB);
1846 molb->nposres_xB = 0;
1848 for(i=0; i<il->nr; i+=2) {
1849 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1850 a = il->iatoms[i+1];
1851 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
1852 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
1853 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
1855 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
1856 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
1857 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
1862 static void set_disres_npair(gmx_mtop_t *mtop)
1869 ip = mtop->ffparams.iparams;
1871 for(mt=0; mt<mtop->nmoltype; mt++) {
1872 il = &mtop->moltype[mt].ilist[F_DISRES];
1876 for(i=0; i<il->nr; i+=3) {
1878 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
1879 ip[a[i]].disres.npair = npair;
1887 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
1895 do_symtab(fio, &(mtop->symtab),bRead);
1897 pr_symtab(debug,0,"symtab",&mtop->symtab);
1899 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
1901 if (file_version >= 57) {
1902 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
1904 gmx_fio_do_int(fio,mtop->nmoltype);
1909 snew(mtop->moltype,mtop->nmoltype);
1910 if (file_version < 57) {
1911 mtop->moltype[0].name = mtop->name;
1914 for(mt=0; mt<mtop->nmoltype; mt++) {
1915 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
1919 if (file_version >= 57) {
1920 gmx_fio_do_int(fio,mtop->nmolblock);
1922 mtop->nmolblock = 1;
1925 snew(mtop->molblock,mtop->nmolblock);
1927 if (file_version >= 57) {
1928 for(mb=0; mb<mtop->nmolblock; mb++) {
1929 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
1931 gmx_fio_do_int(fio,mtop->natoms);
1933 mtop->molblock[0].type = 0;
1934 mtop->molblock[0].nmol = 1;
1935 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
1936 mtop->molblock[0].nposres_xA = 0;
1937 mtop->molblock[0].nposres_xB = 0;
1940 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
1942 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
1944 if (file_version < 57) {
1945 /* Debug statements are inside do_idef */
1946 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
1947 mtop->natoms = mtop->moltype[0].atoms.nr;
1950 if(file_version >= 65)
1952 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
1956 mtop->ffparams.cmap_grid.ngrid = 0;
1957 mtop->ffparams.cmap_grid.grid_spacing = 0.1;
1958 mtop->ffparams.cmap_grid.cmapdata = NULL;
1961 if (file_version >= 57) {
1962 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
1965 if (file_version < 57) {
1966 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
1967 if (bRead && gmx_debug_at) {
1968 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
1970 do_block(fio, &mtop->mols,bRead,file_version);
1971 /* Add the posres coordinates to the molblock */
1972 add_posres_molblock(mtop);
1975 if (file_version >= 57) {
1976 mtop->mols = mtop_mols(mtop);
1979 pr_block(debug,0,"mols",&mtop->mols,TRUE);
1983 if (file_version < 51) {
1984 /* Here used to be the shake blocks */
1985 do_blocka(fio, &dumb,bRead,file_version);
1993 close_symtab(&(mtop->symtab));
1997 /* If TopOnlyOK is TRUE then we can read even future versions
1998 * of tpx files, provided the file_generation hasn't changed.
1999 * If it is FALSE, we need the inputrecord too, and bail out
2000 * if the file is newer than the program.
2002 * The version and generation if the topology (see top of this file)
2003 * are returned in the two last arguments.
2005 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2007 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2008 gmx_bool TopOnlyOK, int *file_version,
2009 int *file_generation)
2018 gmx_fio_checktype(fio);
2019 gmx_fio_setdebug(fio,bDebugMode());
2021 /* NEW! XDR tpb file */
2022 precision = sizeof(real);
2024 gmx_fio_do_string(fio,buf);
2025 if (strncmp(buf,"VERSION",7))
2026 gmx_fatal(FARGS,"Can not read file %s,\n"
2027 " this file is from a Gromacs version which is older than 2.0\n"
2028 " Make a new one with grompp or use a gro or pdb file, if possible",
2029 gmx_fio_getname(fio));
2030 gmx_fio_do_int(fio,precision);
2031 bDouble = (precision == sizeof(double));
2032 if ((precision != sizeof(float)) && !bDouble)
2033 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2034 "instead of %d or %d",
2035 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2036 gmx_fio_setprecision(fio,bDouble);
2037 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2038 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2041 gmx_fio_write_string(fio,GromacsVersion());
2042 bDouble = (precision == sizeof(double));
2043 gmx_fio_setprecision(fio,bDouble);
2044 gmx_fio_do_int(fio,precision);
2046 fgen = tpx_generation;
2049 /* Check versions! */
2050 gmx_fio_do_int(fio,fver);
2053 gmx_fio_do_int(fio,fgen);
2057 if(file_version!=NULL)
2058 *file_version = fver;
2059 if(file_generation!=NULL)
2060 *file_generation = fgen;
2063 if ((fver <= tpx_incompatible_version) ||
2064 ((fver > tpx_version) && !TopOnlyOK) ||
2065 (fgen > tpx_generation))
2066 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2067 gmx_fio_getname(fio),fver,tpx_version);
2069 do_section(fio,eitemHEADER,bRead);
2070 gmx_fio_do_int(fio,tpx->natoms);
2072 gmx_fio_do_int(fio,tpx->ngtc);
2076 gmx_fio_do_int(fio,idum);
2077 gmx_fio_do_real(fio,rdum);
2079 gmx_fio_do_real(fio,tpx->lambda);
2080 gmx_fio_do_int(fio,tpx->bIr);
2081 gmx_fio_do_int(fio,tpx->bTop);
2082 gmx_fio_do_int(fio,tpx->bX);
2083 gmx_fio_do_int(fio,tpx->bV);
2084 gmx_fio_do_int(fio,tpx->bF);
2085 gmx_fio_do_int(fio,tpx->bBox);
2087 if((fgen > tpx_generation)) {
2088 /* This can only happen if TopOnlyOK=TRUE */
2093 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2094 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2095 gmx_bool bXVallocated)
2100 gmx_bool TopOnlyOK,bDum=TRUE;
2101 int file_version,file_generation;
2105 gmx_bool bPeriodicMols;
2108 tpx.natoms = state->natoms;
2109 tpx.ngtc = state->ngtc;
2110 tpx.lambda = state->lambda;
2111 tpx.bIr = (ir != NULL);
2112 tpx.bTop = (mtop != NULL);
2113 tpx.bX = (state->x != NULL);
2114 tpx.bV = (state->v != NULL);
2115 tpx.bF = (f != NULL);
2119 TopOnlyOK = (ir==NULL);
2121 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2125 state->lambda = tpx.lambda;
2126 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2130 init_state(state,0,tpx.ngtc,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2131 state->natoms = tpx.natoms;
2132 state->nalloc = tpx.natoms;
2136 init_state(state,tpx.natoms,tpx.ngtc,0,0); /* nose-hoover chains */
2140 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2142 do_test(fio,tpx.bBox,state->box);
2143 do_section(fio,eitemBOX,bRead);
2145 gmx_fio_ndo_rvec(fio,state->box,DIM);
2146 if (file_version >= 51) {
2147 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2149 /* We initialize box_rel after reading the inputrec */
2150 clear_mat(state->box_rel);
2152 if (file_version >= 28) {
2153 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2154 if (file_version < 56) {
2156 gmx_fio_ndo_rvec(fio,mdum,DIM);
2161 if (state->ngtc > 0 && file_version >= 28) {
2163 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2164 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2165 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2166 snew(dumv,state->ngtc);
2167 if (file_version < 69) {
2168 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2170 /* These used to be the Berendsen tcoupl_lambda's */
2171 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2175 /* Prior to tpx version 26, the inputrec was here.
2176 * I moved it to enable partial forward-compatibility
2177 * for analysis/viewer programs.
2179 if(file_version<26) {
2180 do_test(fio,tpx.bIr,ir);
2181 do_section(fio,eitemIR,bRead);
2184 do_inputrec(fio, ir,bRead,file_version,
2185 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2187 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2190 do_inputrec(fio, &dum_ir,bRead,file_version,
2191 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2193 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2194 done_inputrec(&dum_ir);
2200 do_test(fio,tpx.bTop,mtop);
2201 do_section(fio,eitemTOP,bRead);
2204 do_mtop(fio,mtop,bRead, file_version);
2206 do_mtop(fio,&dum_top,bRead,file_version);
2207 done_mtop(&dum_top,TRUE);
2210 do_test(fio,tpx.bX,state->x);
2211 do_section(fio,eitemX,bRead);
2214 state->flags |= (1<<estX);
2216 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2219 do_test(fio,tpx.bV,state->v);
2220 do_section(fio,eitemV,bRead);
2223 state->flags |= (1<<estV);
2225 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2228 do_test(fio,tpx.bF,f);
2229 do_section(fio,eitemF,bRead);
2230 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2232 /* Starting with tpx version 26, we have the inputrec
2233 * at the end of the file, so we can ignore it
2234 * if the file is never than the software (but still the
2235 * same generation - see comments at the top of this file.
2240 bPeriodicMols = FALSE;
2241 if (file_version >= 26) {
2242 do_test(fio,tpx.bIr,ir);
2243 do_section(fio,eitemIR,bRead);
2245 if (file_version >= 53) {
2246 /* Removed the pbc info from do_inputrec, since we always want it */
2249 bPeriodicMols = ir->bPeriodicMols;
2251 gmx_fio_do_int(fio,ePBC);
2252 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2254 if (file_generation <= tpx_generation && ir) {
2255 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2257 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2258 if (file_version < 51)
2259 set_box_rel(ir,state);
2260 if (file_version < 53) {
2262 bPeriodicMols = ir->bPeriodicMols;
2265 if (bRead && ir && file_version >= 53) {
2266 /* We need to do this after do_inputrec, since that initializes ir */
2268 ir->bPeriodicMols = bPeriodicMols;
2277 if (state->ngtc == 0)
2279 /* Reading old version without tcoupl state data: set it */
2280 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2282 if (tpx.bTop && mtop)
2284 if (file_version < 57)
2286 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2288 ir->eDisre = edrSimple;
2292 ir->eDisre = edrNone;
2295 set_disres_npair(mtop);
2299 if (tpx.bTop && mtop)
2301 gmx_mtop_finalize(mtop);
2304 if (file_version >= 57)
2308 env = getenv("GMX_NOCHARGEGROUPS");
2311 sscanf(env,"%d",&ienv);
2312 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2317 "Will make single atomic charge groups in non-solvent%s\n",
2318 ienv > 1 ? " and solvent" : "");
2319 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2321 fprintf(stderr,"\n");
2329 /************************************************************
2331 * The following routines are the exported ones
2333 ************************************************************/
2335 t_fileio *open_tpx(const char *fn,const char *mode)
2337 return gmx_fio_open(fn,mode);
2340 void close_tpx(t_fileio *fio)
2345 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2346 int *file_version, int *file_generation)
2350 fio = open_tpx(fn,"r");
2351 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2355 void write_tpx_state(const char *fn,
2356 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2360 fio = open_tpx(fn,"w");
2361 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2365 void read_tpx_state(const char *fn,
2366 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2370 fio = open_tpx(fn,"r");
2371 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2375 int read_tpx(const char *fn,
2376 t_inputrec *ir, matrix box,int *natoms,
2377 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2385 fio = open_tpx(fn,"r");
2386 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2388 *natoms = state.natoms;
2390 copy_mat(state.box,box);
2398 int read_tpx_top(const char *fn,
2399 t_inputrec *ir, matrix box,int *natoms,
2400 rvec *x,rvec *v,rvec *f,t_topology *top)
2406 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2408 *top = gmx_mtop_t_to_t_topology(&mtop);
2413 gmx_bool fn2bTPX(const char *file)
2415 switch (fn2ftp(file)) {
2425 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2426 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2429 int natoms,i,version,generation;
2430 gmx_bool bTop,bXNULL;
2432 t_topology *topconv;
2435 bTop = fn2bTPX(infile);
2438 read_tpxheader(infile,&header,TRUE,&version,&generation);
2440 snew(*x,header.natoms);
2442 snew(*v,header.natoms);
2444 *ePBC = read_tpx(infile,NULL,box,&natoms,
2445 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2446 *top = gmx_mtop_t_to_t_topology(mtop);
2448 strcpy(title,*top->name);
2449 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2452 get_stx_coordnum(infile,&natoms);
2453 init_t_atoms(&top->atoms,natoms,FALSE);
2454 bXNULL = (x == NULL);
2458 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2464 aps = gmx_atomprop_init();
2465 for(i=0; (i<natoms); i++)
2466 if (!gmx_atomprop_query(aps,epropMass,
2467 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2468 *top->atoms.atomname[i],
2469 &(top->atoms.atom[i].m))) {
2471 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2472 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2473 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2474 *top->atoms.atomname[i]);
2476 gmx_atomprop_destroy(aps);
2478 top->idef.ntypes=-1;