Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / gmxlib / tpxio.c
blob31091ad40a25f1c6254103f89bb865c9194d3fb3
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
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
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 /* This file is completely threadsafe - keep it that way! */
41 #ifdef GMX_THREADS
42 #include <thread_mpi.h>
43 #endif
46 #include <ctype.h>
47 #include "sysstuff.h"
48 #include "smalloc.h"
49 #include "string2.h"
50 #include "gmx_fatal.h"
51 #include "macros.h"
52 #include "names.h"
53 #include "symtab.h"
54 #include "futil.h"
55 #include "filenm.h"
56 #include "gmxfio.h"
57 #include "topsort.h"
58 #include "tpxio.h"
59 #include "txtdump.h"
60 #include "confio.h"
61 #include "atomprop.h"
62 #include "copyrite.h"
63 #include "vec.h"
64 #include "mtop_util.h"
66 /* This number should be increased whenever the file format changes! */
67 static const int tpx_version = 74;
69 /* This number should only be increased when you edit the TOPOLOGY section
70 * of the tpx format. This way we can maintain forward compatibility too
71 * for all analysis tools and/or external programs that only need to
72 * know the atom/residue names, charges, and bond connectivity.
74 * It first appeared in tpx version 26, when I also moved the inputrecord
75 * to the end of the tpx file, so we can just skip it if we only
76 * want the topology.
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 */
88 typedef struct {
89 int fvnr; /* file version number in which the function type first appeared */
90 int ftype; /* function type */
91 } t_ftupd;
93 /*
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[] = {
99 { 20, F_CUBICBONDS },
100 { 20, F_CONNBONDS },
101 { 20, F_HARMONIC },
102 { 20, F_EQM, },
103 { 22, F_DISRESVIOL },
104 { 22, F_ORIRES },
105 { 22, F_ORIRESDEV },
106 { 26, F_FOURDIHS },
107 { 26, F_PIDIHS },
108 { 26, F_DIHRES },
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 }
114 };*/
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 },
123 { 20, F_CONNBONDS },
124 { 20, F_HARMONIC },
125 { 34, F_FENEBONDS },
126 { 43, F_TABBONDS },
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 },
133 { 43, F_TABANGLES },
134 { 26, F_FOURDIHS },
135 { 26, F_PIDIHS },
136 { 43, F_TABDIHS },
137 { 65, F_CMAP },
138 { 60, F_GB12 },
139 { 61, F_GB13 },
140 { 61, F_GB14 },
141 { 72, F_GBPOL },
142 { 72, F_NPSOLVATION },
143 { 41, F_LJC14_Q },
144 { 41, F_LJC_PAIRS_NB },
145 { 32, F_BHAM_LR },
146 { 32, F_RF_EXCL },
147 { 32, F_COUL_RECIP },
148 { 46, F_DPD },
149 { 30, F_POLARIZATION },
150 { 36, F_THOLE_POL },
151 { 22, F_DISRESVIOL },
152 { 22, F_ORIRES },
153 { 22, F_ORIRESDEV },
154 { 26, F_DIHRES },
155 { 26, F_DIHRESVIOL },
156 { 49, F_VSITE4FDN },
157 { 50, F_VSITEN },
158 { 46, F_COM_PULL },
159 { 20, F_EQM },
160 { 46, F_ECONSERVED },
161 { 69, F_VTEMP },
162 { 66, F_PDISPCORR },
163 { 54, F_DHDL_CON },
165 #define NFTUPD asize(ftupd)
167 /* Needed for backward compatibility */
168 #define MAXNODES 256
170 static void _do_section(t_fileio *fio,int key,gmx_bool bRead,const char *src,
171 int line)
173 char buf[STRLEN];
174 gmx_bool bDbg;
176 if (gmx_fio_getftp(fio) == efTPA) {
177 if (!bRead) {
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);
184 else {
185 if (gmx_fio_getdebug(fio))
186 fprintf(stderr,"Looking for section %s (%s, %d)",
187 itemstr[key],src,line);
189 do {
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,
209 int file_version)
211 gmx_bool bDum=TRUE;
212 int i;
214 gmx_fio_do_int(fio,pgrp->nat);
215 if (bRead)
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);
219 if (bRead)
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);
229 } else {
230 pgrp->kB = pgrp->k;
234 static void do_pull(t_fileio *fio, t_pull *pull,gmx_bool bRead, int file_version)
236 int g;
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);
246 if (bRead)
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)
255 gmx_bool bDum=TRUE;
256 int i;
258 gmx_fio_do_int(fio,rotg->eType);
259 gmx_fio_do_int(fio,rotg->bMassW);
260 gmx_fio_do_int(fio,rotg->nat);
261 if (bRead)
262 snew(rotg->ind,rotg->nat);
263 gmx_fio_ndo_int(fio,rotg->ind,rotg->nat);
264 if (bRead)
265 snew(rotg->x_ref,rotg->nat);
266 gmx_fio_ndo_rvec(fio,rotg->x_ref,rotg->nat);
267 gmx_fio_do_rvec(fio,rotg->vec);
268 gmx_fio_do_rvec(fio,rotg->pivot);
269 gmx_fio_do_real(fio,rotg->rate);
270 gmx_fio_do_real(fio,rotg->k);
271 gmx_fio_do_real(fio,rotg->slab_dist);
272 gmx_fio_do_real(fio,rotg->min_gaussian);
273 gmx_fio_do_real(fio,rotg->eps);
274 gmx_fio_do_int(fio,rotg->eFittype);
277 static void do_rot(t_fileio *fio, t_rot *rot,gmx_bool bRead, int file_version)
279 int g;
281 gmx_fio_do_int(fio,rot->ngrp);
282 gmx_fio_do_int(fio,rot->nstrout);
283 gmx_fio_do_int(fio,rot->nstsout);
284 if (bRead)
285 snew(rot->grp,rot->ngrp);
286 for(g=0; g<rot->ngrp; g++)
287 do_rotgrp(fio, &rot->grp[g],bRead,file_version);
291 static void do_inputrec(t_fileio *fio, t_inputrec *ir,gmx_bool bRead,
292 int file_version, real *fudgeQQ)
294 int i,j,k,*tmp,idum=0;
295 gmx_bool bDum=TRUE;
296 real rdum,bd_temp;
297 rvec vdum;
298 gmx_bool bSimAnn;
299 real zerotemptime,finish_t,init_temp,finish_temp;
301 if (file_version != tpx_version)
303 /* Give a warning about features that are not accessible */
304 fprintf(stderr,"Note: tpx file_version %d, software version %d\n",
305 file_version,tpx_version);
308 if (bRead)
310 init_inputrec(ir);
313 if (file_version == 0)
315 return;
318 /* Basic inputrec stuff */
319 gmx_fio_do_int(fio,ir->eI);
320 if (file_version >= 62) {
321 gmx_fio_do_gmx_large_int(fio, ir->nsteps);
322 } else {
323 gmx_fio_do_int(fio,idum);
324 ir->nsteps = idum;
326 if(file_version > 25) {
327 if (file_version >= 62) {
328 gmx_fio_do_gmx_large_int(fio, ir->init_step);
329 } else {
330 gmx_fio_do_int(fio,idum);
331 ir->init_step = idum;
333 } else {
334 ir->init_step=0;
337 if(file_version >= 58)
338 gmx_fio_do_int(fio,ir->simulation_part);
339 else
340 ir->simulation_part=1;
342 if (file_version >= 67) {
343 gmx_fio_do_int(fio,ir->nstcalcenergy);
344 } else {
345 ir->nstcalcenergy = 1;
347 if (file_version < 53) {
348 /* The pbc info has been moved out of do_inputrec,
349 * since we always want it, also without reading the inputrec.
351 gmx_fio_do_int(fio,ir->ePBC);
352 if ((file_version <= 15) && (ir->ePBC == 2))
353 ir->ePBC = epbcNONE;
354 if (file_version >= 45) {
355 gmx_fio_do_int(fio,ir->bPeriodicMols);
356 } else {
357 if (ir->ePBC == 2) {
358 ir->ePBC = epbcXYZ;
359 ir->bPeriodicMols = TRUE;
360 } else {
361 ir->bPeriodicMols = FALSE;
365 gmx_fio_do_int(fio,ir->ns_type);
366 gmx_fio_do_int(fio,ir->nstlist);
367 gmx_fio_do_int(fio,ir->ndelta);
368 if (file_version < 41) {
369 gmx_fio_do_int(fio,idum);
370 gmx_fio_do_int(fio,idum);
372 if (file_version >= 45)
373 gmx_fio_do_real(fio,ir->rtpi);
374 else
375 ir->rtpi = 0.05;
376 gmx_fio_do_int(fio,ir->nstcomm);
377 if (file_version > 34)
378 gmx_fio_do_int(fio,ir->comm_mode);
379 else if (ir->nstcomm < 0)
380 ir->comm_mode = ecmANGULAR;
381 else
382 ir->comm_mode = ecmLINEAR;
383 ir->nstcomm = abs(ir->nstcomm);
385 if(file_version > 25)
386 gmx_fio_do_int(fio,ir->nstcheckpoint);
387 else
388 ir->nstcheckpoint=0;
390 gmx_fio_do_int(fio,ir->nstcgsteep);
392 if(file_version>=30)
393 gmx_fio_do_int(fio,ir->nbfgscorr);
394 else if (bRead)
395 ir->nbfgscorr = 10;
397 gmx_fio_do_int(fio,ir->nstlog);
398 gmx_fio_do_int(fio,ir->nstxout);
399 gmx_fio_do_int(fio,ir->nstvout);
400 gmx_fio_do_int(fio,ir->nstfout);
401 gmx_fio_do_int(fio,ir->nstenergy);
402 gmx_fio_do_int(fio,ir->nstxtcout);
403 if (file_version >= 59) {
404 gmx_fio_do_double(fio,ir->init_t);
405 gmx_fio_do_double(fio,ir->delta_t);
406 } else {
407 gmx_fio_do_real(fio,rdum);
408 ir->init_t = rdum;
409 gmx_fio_do_real(fio,rdum);
410 ir->delta_t = rdum;
412 gmx_fio_do_real(fio,ir->xtcprec);
413 if (file_version < 19) {
414 gmx_fio_do_int(fio,idum);
415 gmx_fio_do_int(fio,idum);
417 if(file_version < 18)
418 gmx_fio_do_int(fio,idum);
419 gmx_fio_do_real(fio,ir->rlist);
420 if (file_version >= 67) {
421 gmx_fio_do_real(fio,ir->rlistlong);
423 gmx_fio_do_int(fio,ir->coulombtype);
424 if (file_version < 32 && ir->coulombtype == eelRF)
425 ir->coulombtype = eelRF_NEC;
426 gmx_fio_do_real(fio,ir->rcoulomb_switch);
427 gmx_fio_do_real(fio,ir->rcoulomb);
428 gmx_fio_do_int(fio,ir->vdwtype);
429 gmx_fio_do_real(fio,ir->rvdw_switch);
430 gmx_fio_do_real(fio,ir->rvdw);
431 if (file_version < 67) {
432 ir->rlistlong = max_cutoff(ir->rlist,max_cutoff(ir->rvdw,ir->rcoulomb));
434 gmx_fio_do_int(fio,ir->eDispCorr);
435 gmx_fio_do_real(fio,ir->epsilon_r);
436 if (file_version >= 37) {
437 gmx_fio_do_real(fio,ir->epsilon_rf);
438 } else {
439 if (EEL_RF(ir->coulombtype)) {
440 ir->epsilon_rf = ir->epsilon_r;
441 ir->epsilon_r = 1.0;
442 } else {
443 ir->epsilon_rf = 1.0;
446 if (file_version >= 29)
447 gmx_fio_do_real(fio,ir->tabext);
448 else
449 ir->tabext=1.0;
451 if(file_version > 25) {
452 gmx_fio_do_int(fio,ir->gb_algorithm);
453 gmx_fio_do_int(fio,ir->nstgbradii);
454 gmx_fio_do_real(fio,ir->rgbradii);
455 gmx_fio_do_real(fio,ir->gb_saltconc);
456 gmx_fio_do_int(fio,ir->implicit_solvent);
457 } else {
458 ir->gb_algorithm=egbSTILL;
459 ir->nstgbradii=1;
460 ir->rgbradii=1.0;
461 ir->gb_saltconc=0;
462 ir->implicit_solvent=eisNO;
464 if(file_version>=55)
466 gmx_fio_do_real(fio,ir->gb_epsilon_solvent);
467 gmx_fio_do_real(fio,ir->gb_obc_alpha);
468 gmx_fio_do_real(fio,ir->gb_obc_beta);
469 gmx_fio_do_real(fio,ir->gb_obc_gamma);
470 if(file_version>=60)
472 gmx_fio_do_real(fio,ir->gb_dielectric_offset);
473 gmx_fio_do_int(fio,ir->sa_algorithm);
475 else
477 ir->gb_dielectric_offset = 0.009;
478 ir->sa_algorithm = esaAPPROX;
480 gmx_fio_do_real(fio,ir->sa_surface_tension);
482 /* Override sa_surface_tension if it is not changed in the mpd-file */
483 if(ir->sa_surface_tension<0)
485 if(ir->gb_algorithm==egbSTILL)
487 ir->sa_surface_tension = 0.0049 * 100 * CAL2JOULE;
489 else if(ir->gb_algorithm==egbHCT || ir->gb_algorithm==egbOBC)
491 ir->sa_surface_tension = 0.0054 * 100 * CAL2JOULE;
496 else
498 /* Better use sensible values than insane (0.0) ones... */
499 ir->gb_epsilon_solvent = 80;
500 ir->gb_obc_alpha = 1.0;
501 ir->gb_obc_beta = 0.8;
502 ir->gb_obc_gamma = 4.85;
503 ir->sa_surface_tension = 2.092;
507 gmx_fio_do_int(fio,ir->nkx);
508 gmx_fio_do_int(fio,ir->nky);
509 gmx_fio_do_int(fio,ir->nkz);
510 gmx_fio_do_int(fio,ir->pme_order);
511 gmx_fio_do_real(fio,ir->ewald_rtol);
513 if (file_version >=24)
514 gmx_fio_do_int(fio,ir->ewald_geometry);
515 else
516 ir->ewald_geometry=eewg3D;
518 if (file_version <=17) {
519 ir->epsilon_surface=0;
520 if (file_version==17)
521 gmx_fio_do_int(fio,idum);
523 else
524 gmx_fio_do_real(fio,ir->epsilon_surface);
526 gmx_fio_do_gmx_bool(fio,ir->bOptFFT);
528 gmx_fio_do_gmx_bool(fio,ir->bContinuation);
529 gmx_fio_do_int(fio,ir->etc);
530 /* before version 18, ir->etc was a gmx_bool (ir->btc),
531 * but the values 0 and 1 still mean no and
532 * berendsen temperature coupling, respectively.
534 if (file_version >= 71)
536 gmx_fio_do_int(fio,ir->nsttcouple);
538 else
540 ir->nsttcouple = ir->nstcalcenergy;
542 if (file_version <= 15)
544 gmx_fio_do_int(fio,idum);
546 if (file_version <=17)
548 gmx_fio_do_int(fio,ir->epct);
549 if (file_version <= 15)
551 if (ir->epct == 5)
553 ir->epct = epctSURFACETENSION;
555 gmx_fio_do_int(fio,idum);
557 ir->epct -= 1;
558 /* we have removed the NO alternative at the beginning */
559 if(ir->epct==-1)
561 ir->epc=epcNO;
562 ir->epct=epctISOTROPIC;
564 else
566 ir->epc=epcBERENDSEN;
569 else
571 gmx_fio_do_int(fio,ir->epc);
572 gmx_fio_do_int(fio,ir->epct);
574 if (file_version >= 71)
576 gmx_fio_do_int(fio,ir->nstpcouple);
578 else
580 ir->nstpcouple = ir->nstcalcenergy;
582 gmx_fio_do_real(fio,ir->tau_p);
583 if (file_version <= 15) {
584 gmx_fio_do_rvec(fio,vdum);
585 clear_mat(ir->ref_p);
586 for(i=0; i<DIM; i++)
587 ir->ref_p[i][i] = vdum[i];
588 } else {
589 gmx_fio_do_rvec(fio,ir->ref_p[XX]);
590 gmx_fio_do_rvec(fio,ir->ref_p[YY]);
591 gmx_fio_do_rvec(fio,ir->ref_p[ZZ]);
593 if (file_version <= 15) {
594 gmx_fio_do_rvec(fio,vdum);
595 clear_mat(ir->compress);
596 for(i=0; i<DIM; i++)
597 ir->compress[i][i] = vdum[i];
599 else {
600 gmx_fio_do_rvec(fio,ir->compress[XX]);
601 gmx_fio_do_rvec(fio,ir->compress[YY]);
602 gmx_fio_do_rvec(fio,ir->compress[ZZ]);
604 if (file_version >= 47) {
605 gmx_fio_do_int(fio,ir->refcoord_scaling);
606 gmx_fio_do_rvec(fio,ir->posres_com);
607 gmx_fio_do_rvec(fio,ir->posres_comB);
608 } else {
609 ir->refcoord_scaling = erscNO;
610 clear_rvec(ir->posres_com);
611 clear_rvec(ir->posres_comB);
613 if(file_version > 25)
614 gmx_fio_do_int(fio,ir->andersen_seed);
615 else
616 ir->andersen_seed=0;
618 if(file_version < 26) {
619 gmx_fio_do_gmx_bool(fio,bSimAnn);
620 gmx_fio_do_real(fio,zerotemptime);
623 if (file_version < 37)
624 gmx_fio_do_real(fio,rdum);
626 gmx_fio_do_real(fio,ir->shake_tol);
627 if (file_version < 54)
628 gmx_fio_do_real(fio,*fudgeQQ);
629 gmx_fio_do_int(fio,ir->efep);
630 if (file_version <= 14 && ir->efep > efepNO)
631 ir->efep = efepYES;
632 if (file_version >= 59) {
633 gmx_fio_do_double(fio, ir->init_lambda);
634 gmx_fio_do_double(fio, ir->delta_lambda);
635 } else {
636 gmx_fio_do_real(fio,rdum);
637 ir->init_lambda = rdum;
638 gmx_fio_do_real(fio,rdum);
639 ir->delta_lambda = rdum;
641 if (file_version >= 64) {
642 gmx_fio_do_int(fio,ir->n_flambda);
643 if (bRead) {
644 snew(ir->flambda,ir->n_flambda);
646 bDum=gmx_fio_ndo_double(fio,ir->flambda,ir->n_flambda);
647 } else {
648 ir->n_flambda = 0;
649 ir->flambda = NULL;
651 if (file_version >= 13)
652 gmx_fio_do_real(fio,ir->sc_alpha);
653 else
654 ir->sc_alpha = 0;
655 if (file_version >= 38)
656 gmx_fio_do_int(fio,ir->sc_power);
657 else
658 ir->sc_power = 2;
659 if (file_version >= 15)
660 gmx_fio_do_real(fio,ir->sc_sigma);
661 else
662 ir->sc_sigma = 0.3;
663 if (bRead)
665 if (file_version >= 71)
667 ir->sc_sigma_min = ir->sc_sigma;
669 else
671 ir->sc_sigma_min = 0;
674 if (file_version >= 64) {
675 gmx_fio_do_int(fio,ir->nstdhdl);
676 } else {
677 ir->nstdhdl = 1;
680 if (file_version >= 73)
682 gmx_fio_do_int(fio, ir->separate_dhdl_file);
683 gmx_fio_do_int(fio, ir->dhdl_derivatives);
685 else
687 ir->separate_dhdl_file = sepdhdlfileYES;
688 ir->dhdl_derivatives = dhdlderivativesYES;
691 if (file_version >= 71)
693 gmx_fio_do_int(fio,ir->dh_hist_size);
694 gmx_fio_do_double(fio,ir->dh_hist_spacing);
696 else
698 ir->dh_hist_size = 0;
699 ir->dh_hist_spacing = 0.1;
701 if (file_version >= 57) {
702 gmx_fio_do_int(fio,ir->eDisre);
704 gmx_fio_do_int(fio,ir->eDisreWeighting);
705 if (file_version < 22) {
706 if (ir->eDisreWeighting == 0)
707 ir->eDisreWeighting = edrwEqual;
708 else
709 ir->eDisreWeighting = edrwConservative;
711 gmx_fio_do_gmx_bool(fio,ir->bDisreMixed);
712 gmx_fio_do_real(fio,ir->dr_fc);
713 gmx_fio_do_real(fio,ir->dr_tau);
714 gmx_fio_do_int(fio,ir->nstdisreout);
715 if (file_version >= 22) {
716 gmx_fio_do_real(fio,ir->orires_fc);
717 gmx_fio_do_real(fio,ir->orires_tau);
718 gmx_fio_do_int(fio,ir->nstorireout);
719 } else {
720 ir->orires_fc = 0;
721 ir->orires_tau = 0;
722 ir->nstorireout = 0;
724 if(file_version >= 26) {
725 gmx_fio_do_real(fio,ir->dihre_fc);
726 if (file_version < 56) {
727 gmx_fio_do_real(fio,rdum);
728 gmx_fio_do_int(fio,idum);
730 } else {
731 ir->dihre_fc=0;
734 gmx_fio_do_real(fio,ir->em_stepsize);
735 gmx_fio_do_real(fio,ir->em_tol);
736 if (file_version >= 22)
737 gmx_fio_do_gmx_bool(fio,ir->bShakeSOR);
738 else if (bRead)
739 ir->bShakeSOR = TRUE;
740 if (file_version >= 11)
741 gmx_fio_do_int(fio,ir->niter);
742 else if (bRead) {
743 ir->niter = 25;
744 fprintf(stderr,"Note: niter not in run input file, setting it to %d\n",
745 ir->niter);
747 if (file_version >= 21)
748 gmx_fio_do_real(fio,ir->fc_stepsize);
749 else
750 ir->fc_stepsize = 0;
751 gmx_fio_do_int(fio,ir->eConstrAlg);
752 gmx_fio_do_int(fio,ir->nProjOrder);
753 gmx_fio_do_real(fio,ir->LincsWarnAngle);
754 if (file_version <= 14)
755 gmx_fio_do_int(fio,idum);
756 if (file_version >=26)
757 gmx_fio_do_int(fio,ir->nLincsIter);
758 else if (bRead) {
759 ir->nLincsIter = 1;
760 fprintf(stderr,"Note: nLincsIter not in run input file, setting it to %d\n",
761 ir->nLincsIter);
763 if (file_version < 33)
764 gmx_fio_do_real(fio,bd_temp);
765 gmx_fio_do_real(fio,ir->bd_fric);
766 gmx_fio_do_int(fio,ir->ld_seed);
767 if (file_version >= 33) {
768 for(i=0; i<DIM; i++)
769 gmx_fio_do_rvec(fio,ir->deform[i]);
770 } else {
771 for(i=0; i<DIM; i++)
772 clear_rvec(ir->deform[i]);
774 if (file_version >= 14)
775 gmx_fio_do_real(fio,ir->cos_accel);
776 else if (bRead)
777 ir->cos_accel = 0;
778 gmx_fio_do_int(fio,ir->userint1);
779 gmx_fio_do_int(fio,ir->userint2);
780 gmx_fio_do_int(fio,ir->userint3);
781 gmx_fio_do_int(fio,ir->userint4);
782 gmx_fio_do_real(fio,ir->userreal1);
783 gmx_fio_do_real(fio,ir->userreal2);
784 gmx_fio_do_real(fio,ir->userreal3);
785 gmx_fio_do_real(fio,ir->userreal4);
787 /* AdResS stuff */
788 if (file_version >= 68) {
789 gmx_fio_do_int(fio,ir->adress_type);
790 gmx_fio_do_real(fio,ir->adress_const_wf);
791 gmx_fio_do_real(fio,ir->adress_ex_width);
792 gmx_fio_do_real(fio,ir->adress_hy_width);
793 gmx_fio_do_int(fio,ir->adress_icor);
794 gmx_fio_do_int(fio,ir->badress_chempot_dx);
795 gmx_fio_do_int(fio,ir->badress_tf_full_box);
796 gmx_fio_do_int(fio,ir->adress_site);
797 gmx_fio_do_rvec(fio,ir->adress_refs);
798 gmx_fio_do_int(fio,ir->n_adress_tf_grps);
799 gmx_fio_do_real(fio, ir->adress_ex_forcecap);
800 gmx_fio_do_int(fio, ir->n_energy_grps);
801 gmx_fio_do_int(fio,ir->adress_do_hybridpairs);
803 if (bRead)snew(ir->adress_tf_table_index,ir->n_adress_tf_grps);
804 if (ir->n_adress_tf_grps > 0) {
805 bDum=gmx_fio_ndo_int(fio,ir->adress_tf_table_index,ir->n_adress_tf_grps);
807 if (bRead)snew(ir->adress_group_explicit,ir->n_energy_grps);
808 if (ir->n_energy_grps > 0) {
809 bDum=gmx_fio_ndo_int(fio, ir->adress_group_explicit,ir->n_energy_grps);
813 /* pull stuff */
814 if (file_version >= 48) {
815 gmx_fio_do_int(fio,ir->ePull);
816 if (ir->ePull != epullNO) {
817 if (bRead)
818 snew(ir->pull,1);
819 do_pull(fio, ir->pull,bRead,file_version);
821 } else {
822 ir->ePull = epullNO;
825 /* Enforced rotation */
826 if (file_version >= 74) {
827 gmx_fio_do_int(fio,ir->bRot);
828 if (ir->bRot == TRUE) {
829 if (bRead)
830 snew(ir->rot,1);
831 do_rot(fio, ir->rot,bRead,file_version);
833 } else {
834 ir->bRot = FALSE;
837 /* grpopts stuff */
838 gmx_fio_do_int(fio,ir->opts.ngtc);
839 if (file_version >= 69) {
840 gmx_fio_do_int(fio,ir->opts.nhchainlength);
841 } else {
842 ir->opts.nhchainlength = 1;
844 gmx_fio_do_int(fio,ir->opts.ngacc);
845 gmx_fio_do_int(fio,ir->opts.ngfrz);
846 gmx_fio_do_int(fio,ir->opts.ngener);
848 if (bRead) {
849 snew(ir->opts.nrdf, ir->opts.ngtc);
850 snew(ir->opts.ref_t, ir->opts.ngtc);
851 snew(ir->opts.annealing, ir->opts.ngtc);
852 snew(ir->opts.anneal_npoints, ir->opts.ngtc);
853 snew(ir->opts.anneal_time, ir->opts.ngtc);
854 snew(ir->opts.anneal_temp, ir->opts.ngtc);
855 snew(ir->opts.tau_t, ir->opts.ngtc);
856 snew(ir->opts.nFreeze,ir->opts.ngfrz);
857 snew(ir->opts.acc, ir->opts.ngacc);
858 snew(ir->opts.egp_flags,ir->opts.ngener*ir->opts.ngener);
860 if (ir->opts.ngtc > 0) {
861 if (bRead && file_version<13) {
862 snew(tmp,ir->opts.ngtc);
863 bDum=gmx_fio_ndo_int(fio,tmp, ir->opts.ngtc);
864 for(i=0; i<ir->opts.ngtc; i++)
865 ir->opts.nrdf[i] = tmp[i];
866 sfree(tmp);
867 } else {
868 bDum=gmx_fio_ndo_real(fio,ir->opts.nrdf, ir->opts.ngtc);
870 bDum=gmx_fio_ndo_real(fio,ir->opts.ref_t,ir->opts.ngtc);
871 bDum=gmx_fio_ndo_real(fio,ir->opts.tau_t,ir->opts.ngtc);
872 if (file_version<33 && ir->eI==eiBD) {
873 for(i=0; i<ir->opts.ngtc; i++)
874 ir->opts.tau_t[i] = bd_temp;
877 if (ir->opts.ngfrz > 0)
878 bDum=gmx_fio_ndo_ivec(fio,ir->opts.nFreeze,ir->opts.ngfrz);
879 if (ir->opts.ngacc > 0)
880 gmx_fio_ndo_rvec(fio,ir->opts.acc,ir->opts.ngacc);
881 if (file_version >= 12)
882 bDum=gmx_fio_ndo_int(fio,ir->opts.egp_flags,
883 ir->opts.ngener*ir->opts.ngener);
885 if(bRead && file_version < 26) {
886 for(i=0;i<ir->opts.ngtc;i++) {
887 if(bSimAnn) {
888 ir->opts.annealing[i] = eannSINGLE;
889 ir->opts.anneal_npoints[i] = 2;
890 snew(ir->opts.anneal_time[i],2);
891 snew(ir->opts.anneal_temp[i],2);
892 /* calculate the starting/ending temperatures from reft, zerotemptime, and nsteps */
893 finish_t = ir->init_t + ir->nsteps * ir->delta_t;
894 init_temp = ir->opts.ref_t[i]*(1-ir->init_t/zerotemptime);
895 finish_temp = ir->opts.ref_t[i]*(1-finish_t/zerotemptime);
896 ir->opts.anneal_time[i][0] = ir->init_t;
897 ir->opts.anneal_time[i][1] = finish_t;
898 ir->opts.anneal_temp[i][0] = init_temp;
899 ir->opts.anneal_temp[i][1] = finish_temp;
900 } else {
901 ir->opts.annealing[i] = eannNO;
902 ir->opts.anneal_npoints[i] = 0;
905 } else {
906 /* file version 26 or later */
907 /* First read the lists with annealing and npoints for each group */
908 bDum=gmx_fio_ndo_int(fio,ir->opts.annealing,ir->opts.ngtc);
909 bDum=gmx_fio_ndo_int(fio,ir->opts.anneal_npoints,ir->opts.ngtc);
910 for(j=0;j<(ir->opts.ngtc);j++) {
911 k=ir->opts.anneal_npoints[j];
912 if(bRead) {
913 snew(ir->opts.anneal_time[j],k);
914 snew(ir->opts.anneal_temp[j],k);
916 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_time[j],k);
917 bDum=gmx_fio_ndo_real(fio,ir->opts.anneal_temp[j],k);
920 /* Walls */
921 if (file_version >= 45) {
922 gmx_fio_do_int(fio,ir->nwall);
923 gmx_fio_do_int(fio,ir->wall_type);
924 if (file_version >= 50)
925 gmx_fio_do_real(fio,ir->wall_r_linpot);
926 else
927 ir->wall_r_linpot = -1;
928 gmx_fio_do_int(fio,ir->wall_atomtype[0]);
929 gmx_fio_do_int(fio,ir->wall_atomtype[1]);
930 gmx_fio_do_real(fio,ir->wall_density[0]);
931 gmx_fio_do_real(fio,ir->wall_density[1]);
932 gmx_fio_do_real(fio,ir->wall_ewald_zfac);
933 } else {
934 ir->nwall = 0;
935 ir->wall_type = 0;
936 ir->wall_atomtype[0] = -1;
937 ir->wall_atomtype[1] = -1;
938 ir->wall_density[0] = 0;
939 ir->wall_density[1] = 0;
940 ir->wall_ewald_zfac = 3;
942 /* Cosine stuff for electric fields */
943 for(j=0; (j<DIM); j++) {
944 gmx_fio_do_int(fio,ir->ex[j].n);
945 gmx_fio_do_int(fio,ir->et[j].n);
946 if (bRead) {
947 snew(ir->ex[j].a, ir->ex[j].n);
948 snew(ir->ex[j].phi,ir->ex[j].n);
949 snew(ir->et[j].a, ir->et[j].n);
950 snew(ir->et[j].phi,ir->et[j].n);
952 bDum=gmx_fio_ndo_real(fio,ir->ex[j].a, ir->ex[j].n);
953 bDum=gmx_fio_ndo_real(fio,ir->ex[j].phi,ir->ex[j].n);
954 bDum=gmx_fio_ndo_real(fio,ir->et[j].a, ir->et[j].n);
955 bDum=gmx_fio_ndo_real(fio,ir->et[j].phi,ir->et[j].n);
958 /* QMMM stuff */
959 if(file_version>=39){
960 gmx_fio_do_gmx_bool(fio,ir->bQMMM);
961 gmx_fio_do_int(fio,ir->QMMMscheme);
962 gmx_fio_do_real(fio,ir->scalefactor);
963 gmx_fio_do_int(fio,ir->opts.ngQM);
964 if (bRead) {
965 snew(ir->opts.QMmethod, ir->opts.ngQM);
966 snew(ir->opts.QMbasis, ir->opts.ngQM);
967 snew(ir->opts.QMcharge, ir->opts.ngQM);
968 snew(ir->opts.QMmult, ir->opts.ngQM);
969 snew(ir->opts.bSH, ir->opts.ngQM);
970 snew(ir->opts.CASorbitals, ir->opts.ngQM);
971 snew(ir->opts.CASelectrons,ir->opts.ngQM);
972 snew(ir->opts.SAon, ir->opts.ngQM);
973 snew(ir->opts.SAoff, ir->opts.ngQM);
974 snew(ir->opts.SAsteps, ir->opts.ngQM);
975 snew(ir->opts.bOPT, ir->opts.ngQM);
976 snew(ir->opts.bTS, ir->opts.ngQM);
978 if (ir->opts.ngQM > 0) {
979 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmethod,ir->opts.ngQM);
980 bDum=gmx_fio_ndo_int(fio,ir->opts.QMbasis,ir->opts.ngQM);
981 bDum=gmx_fio_ndo_int(fio,ir->opts.QMcharge,ir->opts.ngQM);
982 bDum=gmx_fio_ndo_int(fio,ir->opts.QMmult,ir->opts.ngQM);
983 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bSH,ir->opts.ngQM);
984 bDum=gmx_fio_ndo_int(fio,ir->opts.CASorbitals,ir->opts.ngQM);
985 bDum=gmx_fio_ndo_int(fio,ir->opts.CASelectrons,ir->opts.ngQM);
986 bDum=gmx_fio_ndo_real(fio,ir->opts.SAon,ir->opts.ngQM);
987 bDum=gmx_fio_ndo_real(fio,ir->opts.SAoff,ir->opts.ngQM);
988 bDum=gmx_fio_ndo_int(fio,ir->opts.SAsteps,ir->opts.ngQM);
989 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bOPT,ir->opts.ngQM);
990 bDum=gmx_fio_ndo_gmx_bool(fio,ir->opts.bTS,ir->opts.ngQM);
992 /* end of QMMM stuff */
997 static void do_harm(t_fileio *fio, t_iparams *iparams,gmx_bool bRead)
999 gmx_fio_do_real(fio,iparams->harmonic.rA);
1000 gmx_fio_do_real(fio,iparams->harmonic.krA);
1001 gmx_fio_do_real(fio,iparams->harmonic.rB);
1002 gmx_fio_do_real(fio,iparams->harmonic.krB);
1005 void do_iparams(t_fileio *fio, t_functype ftype,t_iparams *iparams,
1006 gmx_bool bRead, int file_version)
1008 int i;
1009 gmx_bool bDum;
1010 real rdum;
1012 if (!bRead)
1013 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1014 switch (ftype) {
1015 case F_ANGLES:
1016 case F_G96ANGLES:
1017 case F_BONDS:
1018 case F_G96BONDS:
1019 case F_HARMONIC:
1020 case F_IDIHS:
1021 do_harm(fio, iparams,bRead);
1022 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && bRead) {
1023 /* Correct incorrect storage of parameters */
1024 iparams->pdihs.phiB = iparams->pdihs.phiA;
1025 iparams->pdihs.cpB = iparams->pdihs.cpA;
1027 break;
1028 case F_FENEBONDS:
1029 gmx_fio_do_real(fio,iparams->fene.bm);
1030 gmx_fio_do_real(fio,iparams->fene.kb);
1031 break;
1032 case F_RESTRBONDS:
1033 gmx_fio_do_real(fio,iparams->restraint.lowA);
1034 gmx_fio_do_real(fio,iparams->restraint.up1A);
1035 gmx_fio_do_real(fio,iparams->restraint.up2A);
1036 gmx_fio_do_real(fio,iparams->restraint.kA);
1037 gmx_fio_do_real(fio,iparams->restraint.lowB);
1038 gmx_fio_do_real(fio,iparams->restraint.up1B);
1039 gmx_fio_do_real(fio,iparams->restraint.up2B);
1040 gmx_fio_do_real(fio,iparams->restraint.kB);
1041 break;
1042 case F_TABBONDS:
1043 case F_TABBONDSNC:
1044 case F_TABANGLES:
1045 case F_TABDIHS:
1046 gmx_fio_do_real(fio,iparams->tab.kA);
1047 gmx_fio_do_int(fio,iparams->tab.table);
1048 gmx_fio_do_real(fio,iparams->tab.kB);
1049 break;
1050 case F_CROSS_BOND_BONDS:
1051 gmx_fio_do_real(fio,iparams->cross_bb.r1e);
1052 gmx_fio_do_real(fio,iparams->cross_bb.r2e);
1053 gmx_fio_do_real(fio,iparams->cross_bb.krr);
1054 break;
1055 case F_CROSS_BOND_ANGLES:
1056 gmx_fio_do_real(fio,iparams->cross_ba.r1e);
1057 gmx_fio_do_real(fio,iparams->cross_ba.r2e);
1058 gmx_fio_do_real(fio,iparams->cross_ba.r3e);
1059 gmx_fio_do_real(fio,iparams->cross_ba.krt);
1060 break;
1061 case F_UREY_BRADLEY:
1062 gmx_fio_do_real(fio,iparams->u_b.theta);
1063 gmx_fio_do_real(fio,iparams->u_b.ktheta);
1064 gmx_fio_do_real(fio,iparams->u_b.r13);
1065 gmx_fio_do_real(fio,iparams->u_b.kUB);
1066 break;
1067 case F_QUARTIC_ANGLES:
1068 gmx_fio_do_real(fio,iparams->qangle.theta);
1069 bDum=gmx_fio_ndo_real(fio,iparams->qangle.c,5);
1070 break;
1071 case F_BHAM:
1072 gmx_fio_do_real(fio,iparams->bham.a);
1073 gmx_fio_do_real(fio,iparams->bham.b);
1074 gmx_fio_do_real(fio,iparams->bham.c);
1075 break;
1076 case F_MORSE:
1077 gmx_fio_do_real(fio,iparams->morse.b0);
1078 gmx_fio_do_real(fio,iparams->morse.cb);
1079 gmx_fio_do_real(fio,iparams->morse.beta);
1080 break;
1081 case F_CUBICBONDS:
1082 gmx_fio_do_real(fio,iparams->cubic.b0);
1083 gmx_fio_do_real(fio,iparams->cubic.kb);
1084 gmx_fio_do_real(fio,iparams->cubic.kcub);
1085 break;
1086 case F_CONNBONDS:
1087 break;
1088 case F_POLARIZATION:
1089 gmx_fio_do_real(fio,iparams->polarize.alpha);
1090 break;
1091 case F_WATER_POL:
1092 if (file_version < 31)
1093 gmx_fatal(FARGS,"Old tpr files with water_polarization not supported. Make a new.");
1094 gmx_fio_do_real(fio,iparams->wpol.al_x);
1095 gmx_fio_do_real(fio,iparams->wpol.al_y);
1096 gmx_fio_do_real(fio,iparams->wpol.al_z);
1097 gmx_fio_do_real(fio,iparams->wpol.rOH);
1098 gmx_fio_do_real(fio,iparams->wpol.rHH);
1099 gmx_fio_do_real(fio,iparams->wpol.rOD);
1100 break;
1101 case F_THOLE_POL:
1102 gmx_fio_do_real(fio,iparams->thole.a);
1103 gmx_fio_do_real(fio,iparams->thole.alpha1);
1104 gmx_fio_do_real(fio,iparams->thole.alpha2);
1105 gmx_fio_do_real(fio,iparams->thole.rfac);
1106 break;
1107 case F_LJ:
1108 gmx_fio_do_real(fio,iparams->lj.c6);
1109 gmx_fio_do_real(fio,iparams->lj.c12);
1110 break;
1111 case F_LJ14:
1112 gmx_fio_do_real(fio,iparams->lj14.c6A);
1113 gmx_fio_do_real(fio,iparams->lj14.c12A);
1114 gmx_fio_do_real(fio,iparams->lj14.c6B);
1115 gmx_fio_do_real(fio,iparams->lj14.c12B);
1116 break;
1117 case F_LJC14_Q:
1118 gmx_fio_do_real(fio,iparams->ljc14.fqq);
1119 gmx_fio_do_real(fio,iparams->ljc14.qi);
1120 gmx_fio_do_real(fio,iparams->ljc14.qj);
1121 gmx_fio_do_real(fio,iparams->ljc14.c6);
1122 gmx_fio_do_real(fio,iparams->ljc14.c12);
1123 break;
1124 case F_LJC_PAIRS_NB:
1125 gmx_fio_do_real(fio,iparams->ljcnb.qi);
1126 gmx_fio_do_real(fio,iparams->ljcnb.qj);
1127 gmx_fio_do_real(fio,iparams->ljcnb.c6);
1128 gmx_fio_do_real(fio,iparams->ljcnb.c12);
1129 break;
1130 case F_PDIHS:
1131 case F_PIDIHS:
1132 case F_ANGRES:
1133 case F_ANGRESZ:
1134 gmx_fio_do_real(fio,iparams->pdihs.phiA);
1135 gmx_fio_do_real(fio,iparams->pdihs.cpA);
1136 if ((ftype == F_ANGRES || ftype == F_ANGRESZ) && file_version < 42) {
1137 /* Read the incorrectly stored multiplicity */
1138 gmx_fio_do_real(fio,iparams->harmonic.rB);
1139 gmx_fio_do_real(fio,iparams->harmonic.krB);
1140 iparams->pdihs.phiB = iparams->pdihs.phiA;
1141 iparams->pdihs.cpB = iparams->pdihs.cpA;
1142 } else {
1143 gmx_fio_do_real(fio,iparams->pdihs.phiB);
1144 gmx_fio_do_real(fio,iparams->pdihs.cpB);
1145 gmx_fio_do_int(fio,iparams->pdihs.mult);
1147 break;
1148 case F_DISRES:
1149 gmx_fio_do_int(fio,iparams->disres.label);
1150 gmx_fio_do_int(fio,iparams->disres.type);
1151 gmx_fio_do_real(fio,iparams->disres.low);
1152 gmx_fio_do_real(fio,iparams->disres.up1);
1153 gmx_fio_do_real(fio,iparams->disres.up2);
1154 gmx_fio_do_real(fio,iparams->disres.kfac);
1155 break;
1156 case F_ORIRES:
1157 gmx_fio_do_int(fio,iparams->orires.ex);
1158 gmx_fio_do_int(fio,iparams->orires.label);
1159 gmx_fio_do_int(fio,iparams->orires.power);
1160 gmx_fio_do_real(fio,iparams->orires.c);
1161 gmx_fio_do_real(fio,iparams->orires.obs);
1162 gmx_fio_do_real(fio,iparams->orires.kfac);
1163 break;
1164 case F_DIHRES:
1165 gmx_fio_do_int(fio,iparams->dihres.power);
1166 gmx_fio_do_int(fio,iparams->dihres.label);
1167 gmx_fio_do_real(fio,iparams->dihres.phi);
1168 gmx_fio_do_real(fio,iparams->dihres.dphi);
1169 gmx_fio_do_real(fio,iparams->dihres.kfac);
1170 break;
1171 case F_POSRES:
1172 gmx_fio_do_rvec(fio,iparams->posres.pos0A);
1173 gmx_fio_do_rvec(fio,iparams->posres.fcA);
1174 if (bRead && file_version < 27) {
1175 copy_rvec(iparams->posres.pos0A,iparams->posres.pos0B);
1176 copy_rvec(iparams->posres.fcA,iparams->posres.fcB);
1177 } else {
1178 gmx_fio_do_rvec(fio,iparams->posres.pos0B);
1179 gmx_fio_do_rvec(fio,iparams->posres.fcB);
1181 break;
1182 case F_RBDIHS:
1183 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA,NR_RBDIHS);
1184 if(file_version>=25)
1185 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB,NR_RBDIHS);
1186 break;
1187 case F_FOURDIHS:
1188 /* Fourier dihedrals are internally represented
1189 * as Ryckaert-Bellemans since those are faster to compute.
1191 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcA, NR_RBDIHS);
1192 bDum=gmx_fio_ndo_real(fio,iparams->rbdihs.rbcB, NR_RBDIHS);
1193 break;
1194 case F_CONSTR:
1195 case F_CONSTRNC:
1196 gmx_fio_do_real(fio,iparams->constr.dA);
1197 gmx_fio_do_real(fio,iparams->constr.dB);
1198 break;
1199 case F_SETTLE:
1200 gmx_fio_do_real(fio,iparams->settle.doh);
1201 gmx_fio_do_real(fio,iparams->settle.dhh);
1202 break;
1203 case F_VSITE2:
1204 gmx_fio_do_real(fio,iparams->vsite.a);
1205 break;
1206 case F_VSITE3:
1207 case F_VSITE3FD:
1208 case F_VSITE3FAD:
1209 gmx_fio_do_real(fio,iparams->vsite.a);
1210 gmx_fio_do_real(fio,iparams->vsite.b);
1211 break;
1212 case F_VSITE3OUT:
1213 case F_VSITE4FD:
1214 case F_VSITE4FDN:
1215 gmx_fio_do_real(fio,iparams->vsite.a);
1216 gmx_fio_do_real(fio,iparams->vsite.b);
1217 gmx_fio_do_real(fio,iparams->vsite.c);
1218 break;
1219 case F_VSITEN:
1220 gmx_fio_do_int(fio,iparams->vsiten.n);
1221 gmx_fio_do_real(fio,iparams->vsiten.a);
1222 break;
1223 case F_GB12:
1224 case F_GB13:
1225 case F_GB14:
1226 /* We got rid of some parameters in version 68 */
1227 if(bRead && file_version<68)
1229 gmx_fio_do_real(fio,rdum);
1230 gmx_fio_do_real(fio,rdum);
1231 gmx_fio_do_real(fio,rdum);
1232 gmx_fio_do_real(fio,rdum);
1234 gmx_fio_do_real(fio,iparams->gb.sar);
1235 gmx_fio_do_real(fio,iparams->gb.st);
1236 gmx_fio_do_real(fio,iparams->gb.pi);
1237 gmx_fio_do_real(fio,iparams->gb.gbr);
1238 gmx_fio_do_real(fio,iparams->gb.bmlt);
1239 break;
1240 case F_CMAP:
1241 gmx_fio_do_int(fio,iparams->cmap.cmapA);
1242 gmx_fio_do_int(fio,iparams->cmap.cmapB);
1243 break;
1244 default:
1245 gmx_fatal(FARGS,"unknown function type %d (%s) in %s line %d",
1247 ftype,interaction_function[ftype].name,__FILE__,__LINE__);
1249 if (!bRead)
1250 gmx_fio_unset_comment(fio);
1253 static void do_ilist(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,int file_version,
1254 int ftype)
1256 int i,k,idum;
1257 gmx_bool bDum=TRUE;
1259 if (!bRead) {
1260 gmx_fio_set_comment(fio, interaction_function[ftype].name);
1262 if (file_version < 44) {
1263 for(i=0; i<MAXNODES; i++)
1264 gmx_fio_do_int(fio,idum);
1266 gmx_fio_do_int(fio,ilist->nr);
1267 if (bRead)
1268 snew(ilist->iatoms,ilist->nr);
1269 bDum=gmx_fio_ndo_int(fio,ilist->iatoms,ilist->nr);
1270 if (!bRead)
1271 gmx_fio_unset_comment(fio);
1274 static void do_ffparams(t_fileio *fio, gmx_ffparams_t *ffparams,
1275 gmx_bool bRead, int file_version)
1277 int idum,i,j;
1278 gmx_bool bDum=TRUE;
1279 unsigned int k;
1281 gmx_fio_do_int(fio,ffparams->atnr);
1282 if (file_version < 57) {
1283 gmx_fio_do_int(fio,idum);
1285 gmx_fio_do_int(fio,ffparams->ntypes);
1286 if (bRead && debug)
1287 fprintf(debug,"ffparams->atnr = %d, ntypes = %d\n",
1288 ffparams->atnr,ffparams->ntypes);
1289 if (bRead) {
1290 snew(ffparams->functype,ffparams->ntypes);
1291 snew(ffparams->iparams,ffparams->ntypes);
1293 /* Read/write all the function types */
1294 bDum=gmx_fio_ndo_int(fio,ffparams->functype,ffparams->ntypes);
1295 if (bRead && debug)
1296 pr_ivec(debug,0,"functype",ffparams->functype,ffparams->ntypes,TRUE);
1298 if (file_version >= 66) {
1299 gmx_fio_do_double(fio,ffparams->reppow);
1300 } else {
1301 ffparams->reppow = 12.0;
1304 if (file_version >= 57) {
1305 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1308 /* Check whether all these function types are supported by the code.
1309 * In practice the code is backwards compatible, which means that the
1310 * numbering may have to be altered from old numbering to new numbering
1312 for (i=0; (i<ffparams->ntypes); i++) {
1313 if (bRead)
1314 /* Loop over file versions */
1315 for (k=0; (k<NFTUPD); k++)
1316 /* Compare the read file_version to the update table */
1317 if ((file_version < ftupd[k].fvnr) &&
1318 (ffparams->functype[i] >= ftupd[k].ftype)) {
1319 ffparams->functype[i] += 1;
1320 if (debug) {
1321 fprintf(debug,"Incrementing function type %d to %d (due to %s)\n",
1322 i,ffparams->functype[i],
1323 interaction_function[ftupd[k].ftype].longname);
1324 fflush(debug);
1328 do_iparams(fio, ffparams->functype[i],&ffparams->iparams[i],bRead,
1329 file_version);
1330 if (bRead && debug)
1331 pr_iparams(debug,ffparams->functype[i],&ffparams->iparams[i]);
1335 static void do_ilists(t_fileio *fio, t_ilist *ilist,gmx_bool bRead,
1336 int file_version)
1338 int i,j,renum[F_NRE];
1339 gmx_bool bDum=TRUE,bClear;
1340 unsigned int k;
1342 for(j=0; (j<F_NRE); j++) {
1343 bClear = FALSE;
1344 if (bRead)
1345 for (k=0; k<NFTUPD; k++)
1346 if ((file_version < ftupd[k].fvnr) && (j == ftupd[k].ftype))
1347 bClear = TRUE;
1348 if (bClear) {
1349 ilist[j].nr = 0;
1350 ilist[j].iatoms = NULL;
1351 } else {
1352 do_ilist(fio, &ilist[j],bRead,file_version,j);
1355 if (bRead && gmx_debug_at)
1356 pr_ilist(debug,0,interaction_function[j].longname,
1357 functype,&ilist[j],TRUE);
1362 static void do_idef(t_fileio *fio, gmx_ffparams_t *ffparams,gmx_moltype_t *molt,
1363 gmx_bool bRead, int file_version)
1365 do_ffparams(fio, ffparams,bRead,file_version);
1367 if (file_version >= 54) {
1368 gmx_fio_do_real(fio,ffparams->fudgeQQ);
1371 do_ilists(fio, molt->ilist,bRead,file_version);
1374 static void do_block(t_fileio *fio, t_block *block,gmx_bool bRead,int file_version)
1376 int i,idum,dum_nra,*dum_a;
1377 gmx_bool bDum=TRUE;
1379 if (file_version < 44)
1380 for(i=0; i<MAXNODES; i++)
1381 gmx_fio_do_int(fio,idum);
1382 gmx_fio_do_int(fio,block->nr);
1383 if (file_version < 51)
1384 gmx_fio_do_int(fio,dum_nra);
1385 if (bRead) {
1386 block->nalloc_index = block->nr+1;
1387 snew(block->index,block->nalloc_index);
1389 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1391 if (file_version < 51 && dum_nra > 0) {
1392 snew(dum_a,dum_nra);
1393 bDum=gmx_fio_ndo_int(fio,dum_a,dum_nra);
1394 sfree(dum_a);
1398 static void do_blocka(t_fileio *fio, t_blocka *block,gmx_bool bRead,
1399 int file_version)
1401 int i,idum;
1402 gmx_bool bDum=TRUE;
1404 if (file_version < 44)
1405 for(i=0; i<MAXNODES; i++)
1406 gmx_fio_do_int(fio,idum);
1407 gmx_fio_do_int(fio,block->nr);
1408 gmx_fio_do_int(fio,block->nra);
1409 if (bRead) {
1410 block->nalloc_index = block->nr+1;
1411 snew(block->index,block->nalloc_index);
1412 block->nalloc_a = block->nra;
1413 snew(block->a,block->nalloc_a);
1415 bDum=gmx_fio_ndo_int(fio,block->index,block->nr+1);
1416 bDum=gmx_fio_ndo_int(fio,block->a,block->nra);
1419 static void do_atom(t_fileio *fio, t_atom *atom,int ngrp,gmx_bool bRead,
1420 int file_version, gmx_groups_t *groups,int atnr)
1422 int i,myngrp;
1424 gmx_fio_do_real(fio,atom->m);
1425 gmx_fio_do_real(fio,atom->q);
1426 gmx_fio_do_real(fio,atom->mB);
1427 gmx_fio_do_real(fio,atom->qB);
1428 gmx_fio_do_ushort(fio, atom->type);
1429 gmx_fio_do_ushort(fio, atom->typeB);
1430 gmx_fio_do_int(fio,atom->ptype);
1431 gmx_fio_do_int(fio,atom->resind);
1432 if (file_version >= 52)
1433 gmx_fio_do_int(fio,atom->atomnumber);
1434 else if (bRead)
1435 atom->atomnumber = NOTSET;
1436 if (file_version < 23)
1437 myngrp = 8;
1438 else if (file_version < 39)
1439 myngrp = 9;
1440 else
1441 myngrp = ngrp;
1443 if (file_version < 57) {
1444 unsigned char uchar[egcNR];
1445 gmx_fio_ndo_uchar(fio,uchar,myngrp);
1446 for(i=myngrp; (i<ngrp); i++) {
1447 uchar[i] = 0;
1449 /* Copy the old data format to the groups struct */
1450 for(i=0; i<ngrp; i++) {
1451 groups->grpnr[i][atnr] = uchar[i];
1456 static void do_grps(t_fileio *fio, int ngrp,t_grps grps[],gmx_bool bRead,
1457 int file_version)
1459 int i,j,myngrp;
1460 gmx_bool bDum=TRUE;
1462 if (file_version < 23)
1463 myngrp = 8;
1464 else if (file_version < 39)
1465 myngrp = 9;
1466 else
1467 myngrp = ngrp;
1469 for(j=0; (j<ngrp); j++) {
1470 if (j<myngrp) {
1471 gmx_fio_do_int(fio,grps[j].nr);
1472 if (bRead)
1473 snew(grps[j].nm_ind,grps[j].nr);
1474 bDum=gmx_fio_ndo_int(fio,grps[j].nm_ind,grps[j].nr);
1476 else {
1477 grps[j].nr = 1;
1478 snew(grps[j].nm_ind,grps[j].nr);
1483 static void do_symstr(t_fileio *fio, char ***nm,gmx_bool bRead,t_symtab *symtab)
1485 int ls;
1487 if (bRead) {
1488 gmx_fio_do_int(fio,ls);
1489 *nm = get_symtab_handle(symtab,ls);
1491 else {
1492 ls = lookup_symtab(symtab,*nm);
1493 gmx_fio_do_int(fio,ls);
1497 static void do_strstr(t_fileio *fio, int nstr,char ***nm,gmx_bool bRead,
1498 t_symtab *symtab)
1500 int j;
1502 for (j=0; (j<nstr); j++)
1503 do_symstr(fio, &(nm[j]),bRead,symtab);
1506 static void do_resinfo(t_fileio *fio, int n,t_resinfo *ri,gmx_bool bRead,
1507 t_symtab *symtab, int file_version)
1509 int j;
1511 for (j=0; (j<n); j++) {
1512 do_symstr(fio, &(ri[j].name),bRead,symtab);
1513 if (file_version >= 63) {
1514 gmx_fio_do_int(fio,ri[j].nr);
1515 gmx_fio_do_uchar(fio, ri[j].ic);
1516 } else {
1517 ri[j].nr = j + 1;
1518 ri[j].ic = ' ';
1523 static void do_atoms(t_fileio *fio, t_atoms *atoms,gmx_bool bRead,t_symtab *symtab,
1524 int file_version,
1525 gmx_groups_t *groups)
1527 int i;
1529 gmx_fio_do_int(fio,atoms->nr);
1530 gmx_fio_do_int(fio,atoms->nres);
1531 if (file_version < 57) {
1532 gmx_fio_do_int(fio,groups->ngrpname);
1533 for(i=0; i<egcNR; i++) {
1534 groups->ngrpnr[i] = atoms->nr;
1535 snew(groups->grpnr[i],groups->ngrpnr[i]);
1538 if (bRead) {
1539 snew(atoms->atom,atoms->nr);
1540 snew(atoms->atomname,atoms->nr);
1541 snew(atoms->atomtype,atoms->nr);
1542 snew(atoms->atomtypeB,atoms->nr);
1543 snew(atoms->resinfo,atoms->nres);
1544 if (file_version < 57) {
1545 snew(groups->grpname,groups->ngrpname);
1547 atoms->pdbinfo = NULL;
1549 for(i=0; (i<atoms->nr); i++) {
1550 do_atom(fio, &atoms->atom[i],egcNR,bRead, file_version,groups,i);
1552 do_strstr(fio, atoms->nr,atoms->atomname,bRead,symtab);
1553 if (bRead && (file_version <= 20)) {
1554 for(i=0; i<atoms->nr; i++) {
1555 atoms->atomtype[i] = put_symtab(symtab,"?");
1556 atoms->atomtypeB[i] = put_symtab(symtab,"?");
1558 } else {
1559 do_strstr(fio, atoms->nr,atoms->atomtype,bRead,symtab);
1560 do_strstr(fio, atoms->nr,atoms->atomtypeB,bRead,symtab);
1562 do_resinfo(fio, atoms->nres,atoms->resinfo,bRead,symtab,file_version);
1564 if (file_version < 57) {
1565 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1567 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1571 static void do_groups(t_fileio *fio, gmx_groups_t *groups,
1572 gmx_bool bRead,t_symtab *symtab,
1573 int file_version)
1575 int g,n,i;
1576 gmx_bool bDum=TRUE;
1578 do_grps(fio, egcNR,groups->grps,bRead,file_version);
1579 gmx_fio_do_int(fio,groups->ngrpname);
1580 if (bRead) {
1581 snew(groups->grpname,groups->ngrpname);
1583 do_strstr(fio, groups->ngrpname,groups->grpname,bRead,symtab);
1584 for(g=0; g<egcNR; g++) {
1585 gmx_fio_do_int(fio,groups->ngrpnr[g]);
1586 if (groups->ngrpnr[g] == 0) {
1587 if (bRead) {
1588 groups->grpnr[g] = NULL;
1590 } else {
1591 if (bRead) {
1592 snew(groups->grpnr[g],groups->ngrpnr[g]);
1594 bDum=gmx_fio_ndo_uchar(fio, groups->grpnr[g],groups->ngrpnr[g]);
1599 static void do_atomtypes(t_fileio *fio, t_atomtypes *atomtypes,gmx_bool bRead,
1600 t_symtab *symtab,int file_version)
1602 int i,j;
1603 gmx_bool bDum = TRUE;
1605 if (file_version > 25) {
1606 gmx_fio_do_int(fio,atomtypes->nr);
1607 j=atomtypes->nr;
1608 if (bRead) {
1609 snew(atomtypes->radius,j);
1610 snew(atomtypes->vol,j);
1611 snew(atomtypes->surftens,j);
1612 snew(atomtypes->atomnumber,j);
1613 snew(atomtypes->gb_radius,j);
1614 snew(atomtypes->S_hct,j);
1616 bDum=gmx_fio_ndo_real(fio,atomtypes->radius,j);
1617 bDum=gmx_fio_ndo_real(fio,atomtypes->vol,j);
1618 bDum=gmx_fio_ndo_real(fio,atomtypes->surftens,j);
1619 if(file_version >= 40)
1621 bDum=gmx_fio_ndo_int(fio,atomtypes->atomnumber,j);
1623 if(file_version >= 60)
1625 bDum=gmx_fio_ndo_real(fio,atomtypes->gb_radius,j);
1626 bDum=gmx_fio_ndo_real(fio,atomtypes->S_hct,j);
1628 } else {
1629 /* File versions prior to 26 cannot do GBSA,
1630 * so they dont use this structure
1632 atomtypes->nr = 0;
1633 atomtypes->radius = NULL;
1634 atomtypes->vol = NULL;
1635 atomtypes->surftens = NULL;
1636 atomtypes->atomnumber = NULL;
1637 atomtypes->gb_radius = NULL;
1638 atomtypes->S_hct = NULL;
1642 static void do_symtab(t_fileio *fio, t_symtab *symtab,gmx_bool bRead)
1644 int i,nr;
1645 t_symbuf *symbuf;
1646 char buf[STRLEN];
1648 gmx_fio_do_int(fio,symtab->nr);
1649 nr = symtab->nr;
1650 if (bRead) {
1651 snew(symtab->symbuf,1);
1652 symbuf = symtab->symbuf;
1653 symbuf->bufsize = nr;
1654 snew(symbuf->buf,nr);
1655 for (i=0; (i<nr); i++) {
1656 gmx_fio_do_string(fio,buf);
1657 symbuf->buf[i]=strdup(buf);
1660 else {
1661 symbuf = symtab->symbuf;
1662 while (symbuf!=NULL) {
1663 for (i=0; (i<symbuf->bufsize) && (i<nr); i++)
1664 gmx_fio_do_string(fio,symbuf->buf[i]);
1665 nr-=i;
1666 symbuf=symbuf->next;
1668 if (nr != 0)
1669 gmx_fatal(FARGS,"nr of symtab strings left: %d",nr);
1673 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
1675 int i,j,ngrid,gs,nelem;
1677 gmx_fio_do_int(fio,cmap_grid->ngrid);
1678 gmx_fio_do_int(fio,cmap_grid->grid_spacing);
1680 ngrid = cmap_grid->ngrid;
1681 gs = cmap_grid->grid_spacing;
1682 nelem = gs * gs;
1684 if(bRead)
1686 snew(cmap_grid->cmapdata,ngrid);
1688 for(i=0;i<cmap_grid->ngrid;i++)
1690 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1694 for(i=0;i<cmap_grid->ngrid;i++)
1696 for(j=0;j<nelem;j++)
1698 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4]);
1699 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+1]);
1700 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+2]);
1701 gmx_fio_do_real(fio,cmap_grid->cmapdata[i].cmap[j*4+3]);
1707 void tpx_make_chain_identifiers(t_atoms *atoms,t_block *mols)
1709 int m,a,a0,a1,r;
1710 char c,chainid;
1711 int chainnum;
1713 /* We always assign a new chain number, but save the chain id characters
1714 * for larger molecules.
1716 #define CHAIN_MIN_ATOMS 15
1718 chainnum=0;
1719 chainid='A';
1720 for(m=0; m<mols->nr; m++)
1722 a0=mols->index[m];
1723 a1=mols->index[m+1];
1724 if ((a1-a0 >= CHAIN_MIN_ATOMS) && (chainid <= 'Z'))
1726 c=chainid;
1727 chainid++;
1729 else
1731 c=' ';
1733 for(a=a0; a<a1; a++)
1735 atoms->resinfo[atoms->atom[a].resind].chainnum = chainnum;
1736 atoms->resinfo[atoms->atom[a].resind].chainid = c;
1738 chainnum++;
1741 /* Blank out the chain id if there was only one chain */
1742 if (chainid == 'B')
1744 for(r=0; r<atoms->nres; r++)
1746 atoms->resinfo[r].chainid = ' ';
1751 static void do_moltype(t_fileio *fio, gmx_moltype_t *molt,gmx_bool bRead,
1752 t_symtab *symtab, int file_version,
1753 gmx_groups_t *groups)
1755 int i;
1757 if (file_version >= 57) {
1758 do_symstr(fio, &(molt->name),bRead,symtab);
1761 do_atoms(fio, &molt->atoms, bRead, symtab, file_version, groups);
1763 if (bRead && gmx_debug_at) {
1764 pr_atoms(debug,0,"atoms",&molt->atoms,TRUE);
1767 if (file_version >= 57) {
1768 do_ilists(fio, molt->ilist,bRead,file_version);
1770 do_block(fio, &molt->cgs,bRead,file_version);
1771 if (bRead && gmx_debug_at) {
1772 pr_block(debug,0,"cgs",&molt->cgs,TRUE);
1776 /* This used to be in the atoms struct */
1777 do_blocka(fio, &molt->excls, bRead, file_version);
1780 static void do_molblock(t_fileio *fio, gmx_molblock_t *molb,gmx_bool bRead,
1781 int file_version)
1783 int i;
1785 gmx_fio_do_int(fio,molb->type);
1786 gmx_fio_do_int(fio,molb->nmol);
1787 gmx_fio_do_int(fio,molb->natoms_mol);
1788 /* Position restraint coordinates */
1789 gmx_fio_do_int(fio,molb->nposres_xA);
1790 if (molb->nposres_xA > 0) {
1791 if (bRead) {
1792 snew(molb->posres_xA,molb->nposres_xA);
1794 gmx_fio_ndo_rvec(fio,molb->posres_xA,molb->nposres_xA);
1796 gmx_fio_do_int(fio,molb->nposres_xB);
1797 if (molb->nposres_xB > 0) {
1798 if (bRead) {
1799 snew(molb->posres_xB,molb->nposres_xB);
1801 gmx_fio_ndo_rvec(fio,molb->posres_xB,molb->nposres_xB);
1806 static t_block mtop_mols(gmx_mtop_t *mtop)
1808 int mb,m,a,mol;
1809 t_block mols;
1811 mols.nr = 0;
1812 for(mb=0; mb<mtop->nmolblock; mb++) {
1813 mols.nr += mtop->molblock[mb].nmol;
1815 mols.nalloc_index = mols.nr + 1;
1816 snew(mols.index,mols.nalloc_index);
1818 a = 0;
1819 m = 0;
1820 mols.index[m] = a;
1821 for(mb=0; mb<mtop->nmolblock; mb++) {
1822 for(mol=0; mol<mtop->molblock[mb].nmol; mol++) {
1823 a += mtop->molblock[mb].natoms_mol;
1824 m++;
1825 mols.index[m] = a;
1829 return mols;
1832 static void add_posres_molblock(gmx_mtop_t *mtop)
1834 t_ilist *il;
1835 int am,i,mol,a;
1836 gmx_bool bFE;
1837 gmx_molblock_t *molb;
1838 t_iparams *ip;
1840 il = &mtop->moltype[0].ilist[F_POSRES];
1841 if (il->nr == 0) {
1842 return;
1844 am = 0;
1845 bFE = FALSE;
1846 for(i=0; i<il->nr; i+=2) {
1847 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1848 am = max(am,il->iatoms[i+1]);
1849 if (ip->posres.pos0B[XX] != ip->posres.pos0A[XX] ||
1850 ip->posres.pos0B[YY] != ip->posres.pos0A[YY] ||
1851 ip->posres.pos0B[ZZ] != ip->posres.pos0A[ZZ]) {
1852 bFE = TRUE;
1855 /* Make the posres coordinate block end at a molecule end */
1856 mol = 0;
1857 while(am >= mtop->mols.index[mol+1]) {
1858 mol++;
1860 molb = &mtop->molblock[0];
1861 molb->nposres_xA = mtop->mols.index[mol+1];
1862 snew(molb->posres_xA,molb->nposres_xA);
1863 if (bFE) {
1864 molb->nposres_xB = molb->nposres_xA;
1865 snew(molb->posres_xB,molb->nposres_xB);
1866 } else {
1867 molb->nposres_xB = 0;
1869 for(i=0; i<il->nr; i+=2) {
1870 ip = &mtop->ffparams.iparams[il->iatoms[i]];
1871 a = il->iatoms[i+1];
1872 molb->posres_xA[a][XX] = ip->posres.pos0A[XX];
1873 molb->posres_xA[a][YY] = ip->posres.pos0A[YY];
1874 molb->posres_xA[a][ZZ] = ip->posres.pos0A[ZZ];
1875 if (bFE) {
1876 molb->posres_xB[a][XX] = ip->posres.pos0B[XX];
1877 molb->posres_xB[a][YY] = ip->posres.pos0B[YY];
1878 molb->posres_xB[a][ZZ] = ip->posres.pos0B[ZZ];
1883 static void set_disres_npair(gmx_mtop_t *mtop)
1885 int mt,i,npair;
1886 t_iparams *ip;
1887 t_ilist *il;
1888 t_iatom *a;
1890 ip = mtop->ffparams.iparams;
1892 for(mt=0; mt<mtop->nmoltype; mt++) {
1893 il = &mtop->moltype[mt].ilist[F_DISRES];
1894 if (il->nr > 0) {
1895 a = il->iatoms;
1896 npair = 0;
1897 for(i=0; i<il->nr; i+=3) {
1898 npair++;
1899 if (i+3 == il->nr || ip[a[i]].disres.label != ip[a[i+3]].disres.label) {
1900 ip[a[i]].disres.npair = npair;
1901 npair = 0;
1908 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop,gmx_bool bRead,
1909 int file_version)
1911 int mt,mb,i;
1912 t_blocka dumb;
1914 if (bRead)
1915 init_mtop(mtop);
1916 do_symtab(fio, &(mtop->symtab),bRead);
1917 if (bRead && debug)
1918 pr_symtab(debug,0,"symtab",&mtop->symtab);
1920 do_symstr(fio, &(mtop->name),bRead,&(mtop->symtab));
1922 if (file_version >= 57) {
1923 do_ffparams(fio, &mtop->ffparams,bRead,file_version);
1925 gmx_fio_do_int(fio,mtop->nmoltype);
1926 } else {
1927 mtop->nmoltype = 1;
1929 if (bRead) {
1930 snew(mtop->moltype,mtop->nmoltype);
1931 if (file_version < 57) {
1932 mtop->moltype[0].name = mtop->name;
1935 for(mt=0; mt<mtop->nmoltype; mt++) {
1936 do_moltype(fio, &mtop->moltype[mt],bRead,&mtop->symtab,file_version,
1937 &mtop->groups);
1940 if (file_version >= 57) {
1941 gmx_fio_do_int(fio,mtop->nmolblock);
1942 } else {
1943 mtop->nmolblock = 1;
1945 if (bRead) {
1946 snew(mtop->molblock,mtop->nmolblock);
1948 if (file_version >= 57) {
1949 for(mb=0; mb<mtop->nmolblock; mb++) {
1950 do_molblock(fio, &mtop->molblock[mb],bRead,file_version);
1952 gmx_fio_do_int(fio,mtop->natoms);
1953 } else {
1954 mtop->molblock[0].type = 0;
1955 mtop->molblock[0].nmol = 1;
1956 mtop->molblock[0].natoms_mol = mtop->moltype[0].atoms.nr;
1957 mtop->molblock[0].nposres_xA = 0;
1958 mtop->molblock[0].nposres_xB = 0;
1961 do_atomtypes (fio, &(mtop->atomtypes),bRead,&(mtop->symtab), file_version);
1962 if (bRead && debug)
1963 pr_atomtypes(debug,0,"atomtypes",&mtop->atomtypes,TRUE);
1965 if (file_version < 57) {
1966 /* Debug statements are inside do_idef */
1967 do_idef (fio, &mtop->ffparams,&mtop->moltype[0],bRead,file_version);
1968 mtop->natoms = mtop->moltype[0].atoms.nr;
1971 if(file_version >= 65)
1973 do_cmap(fio, &mtop->ffparams.cmap_grid,bRead);
1975 else
1977 mtop->ffparams.cmap_grid.ngrid = 0;
1978 mtop->ffparams.cmap_grid.grid_spacing = 0.1;
1979 mtop->ffparams.cmap_grid.cmapdata = NULL;
1982 if (file_version >= 57) {
1983 do_groups(fio, &mtop->groups,bRead,&(mtop->symtab),file_version);
1986 if (file_version < 57) {
1987 do_block(fio, &mtop->moltype[0].cgs,bRead,file_version);
1988 if (bRead && gmx_debug_at) {
1989 pr_block(debug,0,"cgs",&mtop->moltype[0].cgs,TRUE);
1991 do_block(fio, &mtop->mols,bRead,file_version);
1992 /* Add the posres coordinates to the molblock */
1993 add_posres_molblock(mtop);
1995 if (bRead) {
1996 if (file_version >= 57) {
1997 mtop->mols = mtop_mols(mtop);
1999 if (gmx_debug_at) {
2000 pr_block(debug,0,"mols",&mtop->mols,TRUE);
2004 if (file_version < 51) {
2005 /* Here used to be the shake blocks */
2006 do_blocka(fio, &dumb,bRead,file_version);
2007 if (dumb.nr > 0)
2008 sfree(dumb.index);
2009 if (dumb.nra > 0)
2010 sfree(dumb.a);
2013 if (bRead) {
2014 close_symtab(&(mtop->symtab));
2018 /* If TopOnlyOK is TRUE then we can read even future versions
2019 * of tpx files, provided the file_generation hasn't changed.
2020 * If it is FALSE, we need the inputrecord too, and bail out
2021 * if the file is newer than the program.
2023 * The version and generation if the topology (see top of this file)
2024 * are returned in the two last arguments.
2026 * If possible, we will read the inputrec even when TopOnlyOK is TRUE.
2028 static void do_tpxheader(t_fileio *fio,gmx_bool bRead,t_tpxheader *tpx,
2029 gmx_bool TopOnlyOK, int *file_version,
2030 int *file_generation)
2032 char buf[STRLEN];
2033 gmx_bool bDouble;
2034 int precision;
2035 int fver,fgen;
2036 int idum=0;
2037 real rdum=0;
2039 gmx_fio_checktype(fio);
2040 gmx_fio_setdebug(fio,bDebugMode());
2042 /* NEW! XDR tpb file */
2043 precision = sizeof(real);
2044 if (bRead) {
2045 gmx_fio_do_string(fio,buf);
2046 if (strncmp(buf,"VERSION",7))
2047 gmx_fatal(FARGS,"Can not read file %s,\n"
2048 " this file is from a Gromacs version which is older than 2.0\n"
2049 " Make a new one with grompp or use a gro or pdb file, if possible",
2050 gmx_fio_getname(fio));
2051 gmx_fio_do_int(fio,precision);
2052 bDouble = (precision == sizeof(double));
2053 if ((precision != sizeof(float)) && !bDouble)
2054 gmx_fatal(FARGS,"Unknown precision in file %s: real is %d bytes "
2055 "instead of %d or %d",
2056 gmx_fio_getname(fio),precision,sizeof(float),sizeof(double));
2057 gmx_fio_setprecision(fio,bDouble);
2058 fprintf(stderr,"Reading file %s, %s (%s precision)\n",
2059 gmx_fio_getname(fio),buf,bDouble ? "double" : "single");
2061 else {
2062 gmx_fio_write_string(fio,GromacsVersion());
2063 bDouble = (precision == sizeof(double));
2064 gmx_fio_setprecision(fio,bDouble);
2065 gmx_fio_do_int(fio,precision);
2066 fver = tpx_version;
2067 fgen = tpx_generation;
2070 /* Check versions! */
2071 gmx_fio_do_int(fio,fver);
2073 if(fver>=26)
2074 gmx_fio_do_int(fio,fgen);
2075 else
2076 fgen=0;
2078 if(file_version!=NULL)
2079 *file_version = fver;
2080 if(file_generation!=NULL)
2081 *file_generation = fgen;
2084 if ((fver <= tpx_incompatible_version) ||
2085 ((fver > tpx_version) && !TopOnlyOK) ||
2086 (fgen > tpx_generation))
2087 gmx_fatal(FARGS,"reading tpx file (%s) version %d with version %d program",
2088 gmx_fio_getname(fio),fver,tpx_version);
2090 do_section(fio,eitemHEADER,bRead);
2091 gmx_fio_do_int(fio,tpx->natoms);
2092 if (fver >= 28)
2093 gmx_fio_do_int(fio,tpx->ngtc);
2094 else
2095 tpx->ngtc = 0;
2096 if (fver < 62) {
2097 gmx_fio_do_int(fio,idum);
2098 gmx_fio_do_real(fio,rdum);
2100 gmx_fio_do_real(fio,tpx->lambda);
2101 gmx_fio_do_int(fio,tpx->bIr);
2102 gmx_fio_do_int(fio,tpx->bTop);
2103 gmx_fio_do_int(fio,tpx->bX);
2104 gmx_fio_do_int(fio,tpx->bV);
2105 gmx_fio_do_int(fio,tpx->bF);
2106 gmx_fio_do_int(fio,tpx->bBox);
2108 if((fgen > tpx_generation)) {
2109 /* This can only happen if TopOnlyOK=TRUE */
2110 tpx->bIr=FALSE;
2114 static int do_tpx(t_fileio *fio, gmx_bool bRead,
2115 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop,
2116 gmx_bool bXVallocated)
2118 t_tpxheader tpx;
2119 t_inputrec dum_ir;
2120 gmx_mtop_t dum_top;
2121 gmx_bool TopOnlyOK,bDum=TRUE;
2122 int file_version,file_generation;
2123 int i;
2124 rvec *xptr,*vptr;
2125 int ePBC;
2126 gmx_bool bPeriodicMols;
2128 if (!bRead) {
2129 tpx.natoms = state->natoms;
2130 tpx.ngtc = state->ngtc;
2131 tpx.lambda = state->lambda;
2132 tpx.bIr = (ir != NULL);
2133 tpx.bTop = (mtop != NULL);
2134 tpx.bX = (state->x != NULL);
2135 tpx.bV = (state->v != NULL);
2136 tpx.bF = (f != NULL);
2137 tpx.bBox = TRUE;
2140 TopOnlyOK = (ir==NULL);
2142 do_tpxheader(fio,bRead,&tpx,TopOnlyOK,&file_version,&file_generation);
2144 if (bRead) {
2145 state->flags = 0;
2146 state->lambda = tpx.lambda;
2147 /* The init_state calls initialize the Nose-Hoover xi integrals to zero */
2148 if (bXVallocated) {
2149 xptr = state->x;
2150 vptr = state->v;
2151 init_state(state,0,tpx.ngtc,0,0); /* nose-hoover chains */ /* eventually, need to add nnhpres here? */
2152 state->natoms = tpx.natoms;
2153 state->nalloc = tpx.natoms;
2154 state->x = xptr;
2155 state->v = vptr;
2156 } else {
2157 init_state(state,tpx.natoms,tpx.ngtc,0,0); /* nose-hoover chains */
2161 #define do_test(fio,b,p) if (bRead && (p!=NULL) && !b) gmx_fatal(FARGS,"No %s in %s",#p,gmx_fio_getname(fio))
2163 do_test(fio,tpx.bBox,state->box);
2164 do_section(fio,eitemBOX,bRead);
2165 if (tpx.bBox) {
2166 gmx_fio_ndo_rvec(fio,state->box,DIM);
2167 if (file_version >= 51) {
2168 gmx_fio_ndo_rvec(fio,state->box_rel,DIM);
2169 } else {
2170 /* We initialize box_rel after reading the inputrec */
2171 clear_mat(state->box_rel);
2173 if (file_version >= 28) {
2174 gmx_fio_ndo_rvec(fio,state->boxv,DIM);
2175 if (file_version < 56) {
2176 matrix mdum;
2177 gmx_fio_ndo_rvec(fio,mdum,DIM);
2182 if (state->ngtc > 0 && file_version >= 28) {
2183 real *dumv;
2184 /*ndo_double(state->nosehoover_xi,state->ngtc,bDum);*/
2185 /*ndo_double(state->nosehoover_vxi,state->ngtc,bDum);*/
2186 /*ndo_double(state->therm_integral,state->ngtc,bDum);*/
2187 snew(dumv,state->ngtc);
2188 if (file_version < 69) {
2189 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2191 /* These used to be the Berendsen tcoupl_lambda's */
2192 bDum=gmx_fio_ndo_real(fio,dumv,state->ngtc);
2193 sfree(dumv);
2196 /* Prior to tpx version 26, the inputrec was here.
2197 * I moved it to enable partial forward-compatibility
2198 * for analysis/viewer programs.
2200 if(file_version<26) {
2201 do_test(fio,tpx.bIr,ir);
2202 do_section(fio,eitemIR,bRead);
2203 if (tpx.bIr) {
2204 if (ir) {
2205 do_inputrec(fio, ir,bRead,file_version,
2206 mtop ? &mtop->ffparams.fudgeQQ : NULL);
2207 if (bRead && debug)
2208 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2210 else {
2211 do_inputrec(fio, &dum_ir,bRead,file_version,
2212 mtop ? &mtop->ffparams.fudgeQQ :NULL);
2213 if (bRead && debug)
2214 pr_inputrec(debug,0,"inputrec",&dum_ir,FALSE);
2215 done_inputrec(&dum_ir);
2221 do_test(fio,tpx.bTop,mtop);
2222 do_section(fio,eitemTOP,bRead);
2223 if (tpx.bTop) {
2224 if (mtop) {
2225 do_mtop(fio,mtop,bRead, file_version);
2226 } else {
2227 do_mtop(fio,&dum_top,bRead,file_version);
2228 done_mtop(&dum_top,TRUE);
2231 do_test(fio,tpx.bX,state->x);
2232 do_section(fio,eitemX,bRead);
2233 if (tpx.bX) {
2234 if (bRead) {
2235 state->flags |= (1<<estX);
2237 gmx_fio_ndo_rvec(fio,state->x,state->natoms);
2240 do_test(fio,tpx.bV,state->v);
2241 do_section(fio,eitemV,bRead);
2242 if (tpx.bV) {
2243 if (bRead) {
2244 state->flags |= (1<<estV);
2246 gmx_fio_ndo_rvec(fio,state->v,state->natoms);
2249 do_test(fio,tpx.bF,f);
2250 do_section(fio,eitemF,bRead);
2251 if (tpx.bF) gmx_fio_ndo_rvec(fio,f,state->natoms);
2253 /* Starting with tpx version 26, we have the inputrec
2254 * at the end of the file, so we can ignore it
2255 * if the file is never than the software (but still the
2256 * same generation - see comments at the top of this file.
2260 ePBC = -1;
2261 bPeriodicMols = FALSE;
2262 if (file_version >= 26) {
2263 do_test(fio,tpx.bIr,ir);
2264 do_section(fio,eitemIR,bRead);
2265 if (tpx.bIr) {
2266 if (file_version >= 53) {
2267 /* Removed the pbc info from do_inputrec, since we always want it */
2268 if (!bRead) {
2269 ePBC = ir->ePBC;
2270 bPeriodicMols = ir->bPeriodicMols;
2272 gmx_fio_do_int(fio,ePBC);
2273 gmx_fio_do_gmx_bool(fio,bPeriodicMols);
2275 if (file_generation <= tpx_generation && ir) {
2276 do_inputrec(fio, ir,bRead,file_version,mtop ? &mtop->ffparams.fudgeQQ : NULL);
2277 if (bRead && debug)
2278 pr_inputrec(debug,0,"inputrec",ir,FALSE);
2279 if (file_version < 51)
2280 set_box_rel(ir,state);
2281 if (file_version < 53) {
2282 ePBC = ir->ePBC;
2283 bPeriodicMols = ir->bPeriodicMols;
2286 if (bRead && ir && file_version >= 53) {
2287 /* We need to do this after do_inputrec, since that initializes ir */
2288 ir->ePBC = ePBC;
2289 ir->bPeriodicMols = bPeriodicMols;
2294 if (bRead)
2296 if (tpx.bIr && ir)
2298 if (state->ngtc == 0)
2300 /* Reading old version without tcoupl state data: set it */
2301 init_gtc_state(state,ir->opts.ngtc,0,ir->opts.nhchainlength);
2303 if (tpx.bTop && mtop)
2305 if (file_version < 57)
2307 if (mtop->moltype[0].ilist[F_DISRES].nr > 0)
2309 ir->eDisre = edrSimple;
2311 else
2313 ir->eDisre = edrNone;
2316 set_disres_npair(mtop);
2320 if (tpx.bTop && mtop)
2322 gmx_mtop_finalize(mtop);
2325 if (file_version >= 57)
2327 char *env;
2328 int ienv;
2329 env = getenv("GMX_NOCHARGEGROUPS");
2330 if (env != NULL)
2332 sscanf(env,"%d",&ienv);
2333 fprintf(stderr,"\nFound env.var. GMX_NOCHARGEGROUPS = %d\n",
2334 ienv);
2335 if (ienv > 0)
2337 fprintf(stderr,
2338 "Will make single atomic charge groups in non-solvent%s\n",
2339 ienv > 1 ? " and solvent" : "");
2340 gmx_mtop_make_atomic_charge_groups(mtop,ienv==1);
2342 fprintf(stderr,"\n");
2347 return ePBC;
2350 /************************************************************
2352 * The following routines are the exported ones
2354 ************************************************************/
2356 t_fileio *open_tpx(const char *fn,const char *mode)
2358 return gmx_fio_open(fn,mode);
2361 void close_tpx(t_fileio *fio)
2363 gmx_fio_close(fio);
2366 void read_tpxheader(const char *fn, t_tpxheader *tpx, gmx_bool TopOnlyOK,
2367 int *file_version, int *file_generation)
2369 t_fileio *fio;
2371 fio = open_tpx(fn,"r");
2372 do_tpxheader(fio,TRUE,tpx,TopOnlyOK,file_version,file_generation);
2373 close_tpx(fio);
2376 void write_tpx_state(const char *fn,
2377 t_inputrec *ir,t_state *state,gmx_mtop_t *mtop)
2379 t_fileio *fio;
2381 fio = open_tpx(fn,"w");
2382 do_tpx(fio,FALSE,ir,state,NULL,mtop,FALSE);
2383 close_tpx(fio);
2386 void read_tpx_state(const char *fn,
2387 t_inputrec *ir,t_state *state,rvec *f,gmx_mtop_t *mtop)
2389 t_fileio *fio;
2391 fio = open_tpx(fn,"r");
2392 do_tpx(fio,TRUE,ir,state,f,mtop,FALSE);
2393 close_tpx(fio);
2396 int read_tpx(const char *fn,
2397 t_inputrec *ir, matrix box,int *natoms,
2398 rvec *x,rvec *v,rvec *f,gmx_mtop_t *mtop)
2400 t_fileio *fio;
2401 t_state state;
2402 int ePBC;
2404 state.x = x;
2405 state.v = v;
2406 fio = open_tpx(fn,"r");
2407 ePBC = do_tpx(fio,TRUE,ir,&state,f,mtop,TRUE);
2408 close_tpx(fio);
2409 *natoms = state.natoms;
2410 if (box)
2411 copy_mat(state.box,box);
2412 state.x = NULL;
2413 state.v = NULL;
2414 done_state(&state);
2416 return ePBC;
2419 int read_tpx_top(const char *fn,
2420 t_inputrec *ir, matrix box,int *natoms,
2421 rvec *x,rvec *v,rvec *f,t_topology *top)
2423 gmx_mtop_t mtop;
2424 t_topology *ltop;
2425 int ePBC;
2427 ePBC = read_tpx(fn,ir,box,natoms,x,v,f,&mtop);
2429 *top = gmx_mtop_t_to_t_topology(&mtop);
2431 return ePBC;
2434 gmx_bool fn2bTPX(const char *file)
2436 switch (fn2ftp(file)) {
2437 case efTPR:
2438 case efTPB:
2439 case efTPA:
2440 return TRUE;
2441 default:
2442 return FALSE;
2446 gmx_bool read_tps_conf(const char *infile,char *title,t_topology *top,int *ePBC,
2447 rvec **x,rvec **v,matrix box,gmx_bool bMass)
2449 t_tpxheader header;
2450 int natoms,i,version,generation;
2451 gmx_bool bTop,bXNULL;
2452 gmx_mtop_t *mtop;
2453 t_topology *topconv;
2454 gmx_atomprop_t aps;
2456 bTop = fn2bTPX(infile);
2457 *ePBC = -1;
2458 if (bTop) {
2459 read_tpxheader(infile,&header,TRUE,&version,&generation);
2460 if (x)
2461 snew(*x,header.natoms);
2462 if (v)
2463 snew(*v,header.natoms);
2464 snew(mtop,1);
2465 *ePBC = read_tpx(infile,NULL,box,&natoms,
2466 (x==NULL) ? NULL : *x,(v==NULL) ? NULL : *v,NULL,mtop);
2467 *top = gmx_mtop_t_to_t_topology(mtop);
2468 sfree(mtop);
2469 strcpy(title,*top->name);
2470 tpx_make_chain_identifiers(&top->atoms,&top->mols);
2472 else {
2473 get_stx_coordnum(infile,&natoms);
2474 init_t_atoms(&top->atoms,natoms,FALSE);
2475 bXNULL = (x == NULL);
2476 snew(*x,natoms);
2477 if (v)
2478 snew(*v,natoms);
2479 read_stx_conf(infile,title,&top->atoms,*x,(v==NULL) ? NULL : *v,ePBC,box);
2480 if (bXNULL) {
2481 sfree(*x);
2482 x = NULL;
2484 if (bMass) {
2485 aps = gmx_atomprop_init();
2486 for(i=0; (i<natoms); i++)
2487 if (!gmx_atomprop_query(aps,epropMass,
2488 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2489 *top->atoms.atomname[i],
2490 &(top->atoms.atom[i].m))) {
2491 if (debug)
2492 fprintf(debug,"Can not find mass for atom %s %d %s, setting to 1\n",
2493 *top->atoms.resinfo[top->atoms.atom[i].resind].name,
2494 top->atoms.resinfo[top->atoms.atom[i].resind].nr,
2495 *top->atoms.atomname[i]);
2497 gmx_atomprop_destroy(aps);
2499 top->idef.ntypes=-1;
2502 return bTop;